Mesh Oriented datABase  (version 5.5.0)
An array-based unstructured mesh library
moab::CN Class Reference

Canonical numbering data and functions This class represents canonical ordering of finite-element meshes. Elements in the finite element "zoo" are represented. Canonical numbering denotes the vertex, edge, and face numbers making up each kind of element, and the vertex numbers defining those entities. Functions for evaluating adjacencies and other things based on vertex numbering are also provided. By default, this class defines a zero-based numbering system. For a complete description of this class, see the document "MOAB Canonical Numbering Conventions", Timothy J. Tautges, Sandia National Laboratories Report #SAND2004-xxxx. More...

#include <CN.hpp>

+ Collaboration diagram for moab::CN:

Classes

struct  ConnMap
 
struct  UpConnMap
 

Public Types

enum  { MAX_NODES_PER_ELEMENT = 27 }
 
enum  { MID_EDGE_BIT = 1 << 1 , MID_FACE_BIT = 1 << 2 , MID_REGION_BIT = 1 << 3 }
 
enum  { INTERSECT = 0 , UNION }
 enum used to specify operation type More...
 

Static Public Member Functions

static DimensionPair getDimPair (int entity_type)
 Get the dimension pair corresponding to a dimension. More...
 
static short int GetBasis ()
 get the basis of the numbering system More...
 
static void SetBasis (const int in_basis)
 set the basis of the numbering system More...
 
static const char * EntityTypeName (const EntityType this_type)
 return the string type name for this type More...
 
static EntityType EntityTypeFromName (const char *name)
 given a name, find the corresponding entity type More...
 
static short int Dimension (const EntityType t)
 return the topological entity dimension More...
 
static short int VerticesPerEntity (const EntityType t)
 return the number of (corner) vertices contained in the specified type. More...
 
static short int NumSubEntities (const EntityType t, const int d)
 return the number of sub-entities bounding the entity. More...
 
static EntityType SubEntityType (const EntityType this_type, const int sub_dimension, const int index)
 return the type of a particular sub-entity. More...
 
static void SubEntityVertexIndices (const EntityType this_type, const int sub_dimension, const int sub_index, int sub_entity_conn[])
 return the vertex indices of the specified sub-entity. More...
 
static const short * SubEntityVertexIndices (const EntityType this_type, const int sub_dimension, const int sub_index, EntityType &sub_type, int &num_sub_ent_vertices)
 return the vertex indices of the specified sub-entity. More...
 
static void SubEntityNodeIndices (const EntityType this_topo, const int num_nodes, const int sub_dimension, const int sub_index, EntityType &sub_entity_topo, int &num_sub_entity_nodes, int sub_entity_conn[])
 return the node indices of the specified sub-entity. More...
 
static void SubEntityConn (const void *parent_conn, const EntityType parent_type, const int sub_dimension, const int sub_index, void *sub_entity_conn, int &num_sub_vertices)
 return the vertices of the specified sub entity More...
 
static short int AdjacentSubEntities (const EntityType this_type, const int *source_indices, const int num_source_indices, const int source_dim, const int target_dim, std::vector< int > &index_list, const int operation_type=CN::INTERSECT)
 For a specified set of sides of given dimension, return the intersection or union of all sides of specified target dimension adjacent to those sides. More...
 
static short int SideNumber (const EntityType parent_type, const int *parent_conn, const int *child_conn, const int child_num_verts, const int child_dim, int &side_number, int &sense, int &offset)
 return the side index represented in the input sub-entity connectivity in the input parent entity connectivity array. More...
 
static short int SideNumber (const EntityType parent_type, const unsigned int *parent_conn, const unsigned int *child_conn, const int child_num_verts, const int child_dim, int &side_number, int &sense, int &offset)
 
static short int SideNumber (const EntityType parent_type, const long *parent_conn, const long *child_conn, const int child_num_verts, const int child_dim, int &side_number, int &sense, int &offset)
 
static short int SideNumber (const EntityType parent_type, const unsigned long *parent_conn, const unsigned long *child_conn, const int child_num_verts, const int child_dim, int &side_number, int &sense, int &offset)
 
static short int SideNumber (const EntityType parent_type, const unsigned long long *parent_conn, const unsigned long long *child_conn, const int child_num_verts, const int child_dim, int &side_number, int &sense, int &offset)
 
static short int SideNumber (const EntityType parent_type, void *const *parent_conn, void *const *child_conn, const int child_num_verts, const int child_dim, int &side_number, int &sense, int &offset)
 
static short int SideNumber (const EntityType parent_type, const int *child_conn_indices, const int child_num_verts, const int child_dim, int &side_number, int &sense, int &offset)
 return the side index represented in the input sub-entity connectivity More...
 
static short int OppositeSide (const EntityType parent_type, const int child_index, const int child_dim, int &opposite_index, int &opposite_dim)
 return the dimension and index of the opposite side, given parent entity type and child dimension and index. This function is only defined for certain types of parent/child types: (Parent, Child dim->Opposite dim): (Tri, 1->0), (Tri, 0->1), (Quad, 1->1), (Quad, 0->0), (Tet, 2->0), (Tet, 1->1), (Tet, 0->2), (Hex, 2->2), (Hex, 1->1)(diagonally across element), (Hex, 0->0) (diagonally across element) All other parent types and child dimensions return an error. More...
 
static bool ConnectivityMatch (const int *conn1, const int *conn2, const int num_vertices, int &direct, int &offset)
 given two connectivity arrays, determine whether or not they represent the same entity. More...
 
static bool ConnectivityMatch (const unsigned int *conn1, const unsigned int *conn2, const int num_vertices, int &direct, int &offset)
 
static bool ConnectivityMatch (const long *conn1, const long *conn2, const int num_vertices, int &direct, int &offset)
 
static bool ConnectivityMatch (const unsigned long *conn1, const unsigned long *conn2, const int num_vertices, int &direct, int &offset)
 
static bool ConnectivityMatch (const unsigned long long *conn1, const unsigned long long *conn2, const int num_vertices, int &direct, int &offset)
 
static bool ConnectivityMatch (void *const *conn1, void *const *conn2, const int num_vertices, int &direct, int &offset)
 
static void setPermutation (const EntityType t, const int dim, short int *pvec, const int num_entries, const bool is_reverse=false)
 Set permutation or reverse permutation vector Forward permutation is from CN's numbering into application's ordering; that is, if i is CN's index, pvec[i] is application's index. This function stores the permutation vector for this type and facet dimension, which then is used in calls to permuteThis or revPermuteThis. More...
 
static void resetPermutation (const EntityType t, const int dim)
 Reset permutation or reverse permutation vector. More...
 
static int permuteThis (const EntityType t, const int dim, int *pvec, const int indices_per_ent, const int num_entries)
 Permute a handle array according to permutation vector set with setPermute; permutation is done in-place. More...
 
static int permuteThis (const EntityType t, const int dim, unsigned int *pvec, const int indices_per_ent, const int num_entries)
 
static int permuteThis (const EntityType t, const int dim, long *pvec, const int indices_per_ent, const int num_entries)
 
static int permuteThis (const EntityType t, const int dim, void **pvec, const int indices_per_ent, const int num_entries)
 
static int revPermuteThis (const EntityType t, const int dim, int *pvec, const int indices_per_ent, const int num_entries)
 Reverse permute a handle array according to reverse permutation vector set with setPermute; reverse permutation is done in-place. More...
 
static int revPermuteThis (const EntityType t, const int dim, unsigned int *pvec, const int indices_per_ent, const int num_entries)
 
static int revPermuteThis (const EntityType t, const int dim, long *pvec, const int indices_per_ent, const int num_entries)
 
static int revPermuteThis (const EntityType t, const int dim, void **pvec, const int indices_per_ent, const int num_entries)
 
static bool HasMidEdgeNodes (const EntityType this_type, const int num_verts)
 true if entities of a given type and number of nodes indicates mid edge nodes are present. More...
 
static bool HasMidFaceNodes (const EntityType this_type, const int num_verts)
 true if entities of a given type and number of nodes indicates mid face nodes are present. More...
 
static bool HasMidRegionNodes (const EntityType this_type, const int num_verts)
 true if entities of a given type and number of nodes indicates mid region nodes are present. More...
 
static void HasMidNodes (const EntityType this_type, const int num_verts, int mid_nodes[4])
 true if entities of a given type and number of nodes indicates mid edge/face/region nodes are present. More...
 
static int HasMidNodes (const EntityType this_type, const int num_verts)
 Same as above, except returns a single integer with the bits, from least significant to most significant set to one if the corresponding mid nodes on sub entities of the least dimension (0) to the highest dimension (3) are present in the elment type. More...
 
static void HONodeParent (EntityType elem_type, int num_nodes, int ho_node_index, int &parent_dim, int &parent_index)
 given data about an element and a vertex in that element, return the dimension and index of the sub-entity that the vertex resolves. If it does not resolve a sub-entity, either because it's a corner node or it's not in the element, -1 is returned in both return values. More...
 
static short int HONodeIndex (const EntityType this_type, const int num_verts, const int subfacet_dim, const int subfacet_index)
 for an entity of this type with num_verts vertices, and a specified subfacet (dimension and index), return the index of the higher order node for that entity in this entity's connectivity array More...
 

Static Public Attributes

static MOAB_EXPORT const ConnMap mConnectivityMap [MBMAXTYPE][3]
 
static const UpConnMap mUpConnMap [MBMAXTYPE][4][4]
 
static MOAB_EXPORT const unsigned char midNodesPerType [MBMAXTYPE][MAX_NODES_PER_ELEMENT+1]
 
static short int permuteVec [MBMAXTYPE][3][MAX_SUB_ENTITIES+1]
 Permutation and reverse permutation vectors. More...
 
static short int revPermuteVec [MBMAXTYPE][3][MAX_SUB_ENTITIES+1]
 
static MOAB_EXPORT const DimensionPair TypeDimensionMap []
 this const vector defines the starting and ending EntityType for each dimension, e.g. TypeDimensionMap[2] returns a pair of EntityTypes bounding dimension 2. More...
 

Private Member Functions

 CN ()
 declare private constructor, since we don't want to create any of these More...
 

Static Private Member Functions

static void SwitchBasis (const int old_basis, const int new_basis)
 switch the basis More...
 

Static Private Attributes

static MOAB_EXPORT const char * entityTypeNames []
 entity names More...
 
static MOAB_EXPORT short int numberBasis = 0
 the basis of the numbering system (normally 0 or 1, 0 by default) More...
 
static MOAB_EXPORT short increasingInts []
 

Detailed Description

Canonical numbering data and functions This class represents canonical ordering of finite-element meshes. Elements in the finite element "zoo" are represented. Canonical numbering denotes the vertex, edge, and face numbers making up each kind of element, and the vertex numbers defining those entities. Functions for evaluating adjacencies and other things based on vertex numbering are also provided. By default, this class defines a zero-based numbering system. For a complete description of this class, see the document "MOAB Canonical Numbering Conventions", Timothy J. Tautges, Sandia National Laboratories Report #SAND2004-xxxx.

MOAB, a Mesh-Oriented datABase, is a software component for creating, storing and accessing finite element mesh data.

Copyright 2004 Sandia Corporation. Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains certain rights in this software.

This library is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation; either version 2.1 of the License, or (at your option) any later version.

Author
Tim Tautges
Date
April 2004

Definition at line 56 of file CN.hpp.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
MID_EDGE_BIT 
MID_FACE_BIT 
MID_REGION_BIT 

Definition at line 78 of file CN.hpp.

79  {
80  MID_EDGE_BIT = 1 << 1,
81  MID_FACE_BIT = 1 << 2,
82  MID_REGION_BIT = 1 << 3
83  };

◆ anonymous enum

anonymous enum

enum used to specify operation type

Enumerator
INTERSECT 
UNION 

Definition at line 86 of file CN.hpp.

87  {
88  INTERSECT = 0,
89  UNION
90  };

◆ anonymous enum

anonymous enum
Enumerator
MAX_NODES_PER_ELEMENT 

Definition at line 74 of file CN.hpp.

75  {
77  };

Constructor & Destructor Documentation

◆ CN()

moab::CN::CN ( )
private

declare private constructor, since we don't want to create any of these

Member Function Documentation

◆ AdjacentSubEntities()

short int moab::CN::AdjacentSubEntities ( const EntityType  this_type,
const int *  source_indices,
const int  num_source_indices,
const int  source_dim,
const int  target_dim,
std::vector< int > &  index_list,
const int  operation_type = CN::INTERSECT 
)
static

For a specified set of sides of given dimension, return the intersection or union of all sides of specified target dimension adjacent to those sides.

given an entity and a target dimension & side number, get that entity

Parameters
this_typeType of entity for which sub-entity connectivity is being queried
source_indicesIndices of sides being queried
num_source_indicesNumber of entries in source_indices
source_dimDimension of source entity
target_dimDimension of target entity
index_listIndices of target entities (returned)
operation_typeSpecify either CN::INTERSECT or CN::UNION to get intersection or union of target entity lists over source entities

Definition at line 138 of file CN.cpp.

145 {
146  // first get all the vertex indices
147  std::vector< int > tmp_indices;
148  const int* it1 = source_indices;
149 
150  assert( source_dim <= 3 && target_dim >= 0 && target_dim <= 3 &&
151  // make sure we're not stepping off the end of the array;
152  ( ( source_dim > 0 && *it1 < mConnectivityMap[this_type][source_dim - 1].num_sub_elements ) ||
153  ( source_dim == 0 &&
154  *it1 < mConnectivityMap[this_type][Dimension( this_type ) - 1].num_corners_per_sub_element[0] ) ) &&
155  *it1 >= 0 );
156 
157 #define MUC CN::mUpConnMap[this_type][source_dim][target_dim]
158 
159  // if we're looking for the vertices of a single side, return them in
160  // the canonical ordering; otherwise, return them in sorted order
161  if( num_source_indices == 1 && 0 == target_dim && source_dim != target_dim )
162  {
163 
164  // element of mConnectivityMap should be for this type and for one
165  // less than source_dim, which should give the connectivity of that sub element
166  const ConnMap& cm = mConnectivityMap[this_type][source_dim - 1];
167  std::copy( cm.conn[source_indices[0]],
168  cm.conn[source_indices[0]] + cm.num_corners_per_sub_element[source_indices[0]],
169  std::back_inserter( index_list ) );
170  return 0;
171  }
172 
173  // now go through source indices, folding adjacencies into target list
174  for( it1 = source_indices; it1 != source_indices + num_source_indices; it1++ )
175  {
176  // *it1 is the side index
177  // at start of iteration, index_list has the target list
178 
179  // if a union, or first iteration and index list was empty, copy the list
180  if( operation_type == CN::UNION || ( it1 == source_indices && index_list.empty() ) )
181  {
182  std::copy( MUC.targets_per_source_element[*it1],
183  MUC.targets_per_source_element[*it1] + MUC.num_targets_per_source_element[*it1],
184  std::back_inserter( index_list ) );
185  }
186  else
187  {
188  // else we're intersecting, and have a non-empty list; intersect with this target list
189  tmp_indices.clear();
190  for( int i = MUC.num_targets_per_source_element[*it1] - 1; i >= 0; i-- )
191  if( std::find( index_list.begin(), index_list.end(), MUC.targets_per_source_element[*it1][i] ) !=
192  index_list.end() )
193  tmp_indices.push_back( MUC.targets_per_source_element[*it1][i] );
194  // std::set_intersection(MUC.targets_per_source_element[*it1],
195  // MUC.targets_per_source_element[*it1]+
196  // MUC.num_targets_per_source_element[*it1],
197  // index_list.begin(), index_list.end(),
198  // std::back_inserter(tmp_indices));
199  index_list.swap( tmp_indices );
200 
201  // if we're at this point and the list is empty, the intersection will be NULL;
202  // return if so
203  if( index_list.empty() ) return 0;
204  }
205  }
206 
207  if( operation_type == CN::UNION && num_source_indices != 1 )
208  {
209  // need to sort then unique the list
210  std::sort( index_list.begin(), index_list.end() );
211  index_list.erase( std::unique( index_list.begin(), index_list.end() ), index_list.end() );
212  }
213 
214  return 0;
215 }

References moab::CN::ConnMap::conn, Dimension(), mConnectivityMap, MUC, moab::CN::ConnMap::num_corners_per_sub_element, and UNION.

Referenced by moab::AEntityFactory::get_down_adjacency_elements(), and moab::Core::side_element().

◆ ConnectivityMatch() [1/6]

bool moab::CN::ConnectivityMatch ( const int *  conn1,
const int *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
)
static

given two connectivity arrays, determine whether or not they represent the same entity.

Parameters
conn1Connectivity array of first entity
conn2Connectivity array of second entity
num_verticesNumber of entries in conn1 and conn2
directIf positive, entities have the same sense (returned)
offsetOffset of conn2's first vertex in conn1
Returns
bool Returns true if conn1 and conn2 match

Definition at line 540 of file CN.cpp.

541 {
542  return connectivity_match< int >( conn1_i, conn2_i, num_vertices, direct, offset );
543 }

Referenced by moab::HalfFacetRep::get_down_adjacencies_face_3d(), moab::Core::merge_entities(), moab::Core::side_number(), and SideNumber().

◆ ConnectivityMatch() [2/6]

bool moab::CN::ConnectivityMatch ( const long *  conn1,
const long *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
)
static

Definition at line 554 of file CN.cpp.

555 {
556  return connectivity_match< long >( conn1_i, conn2_i, num_vertices, direct, offset );
557 }

◆ ConnectivityMatch() [3/6]

bool moab::CN::ConnectivityMatch ( const unsigned int *  conn1,
const unsigned int *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
)
static

Definition at line 545 of file CN.cpp.

550 {
551  return connectivity_match< unsigned int >( conn1_i, conn2_i, num_vertices, direct, offset );
552 }

◆ ConnectivityMatch() [4/6]

bool moab::CN::ConnectivityMatch ( const unsigned long *  conn1,
const unsigned long *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
)
static

Definition at line 559 of file CN.cpp.

564 {
565  return connectivity_match< unsigned long >( conn1_i, conn2_i, num_vertices, direct, offset );
566 }

◆ ConnectivityMatch() [5/6]

bool moab::CN::ConnectivityMatch ( const unsigned long long *  conn1,
const unsigned long long *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
)
static

Definition at line 568 of file CN.cpp.

573 {
574  return connectivity_match< unsigned long long >( conn1_i, conn2_i, num_vertices, direct, offset );
575 }

◆ ConnectivityMatch() [6/6]

bool moab::CN::ConnectivityMatch ( void *const *  conn1,
void *const *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
)
static

Definition at line 577 of file CN.cpp.

582 {
583  return connectivity_match< void* >( conn1_i, conn2_i, num_vertices, direct, offset );
584 }

◆ Dimension()

short int moab::CN::Dimension ( const EntityType  t)
static

return the topological entity dimension

Definition at line 1066 of file CN.cpp.

1067 {
1068  return mConnectivityMap[t][0].topo_dimension;
1069 }

References mConnectivityMap, t, and moab::CN::ConnMap::topo_dimension.

Referenced by AdjacentSubEntities(), moab::Range::all_of_dimension(), moab::SweptElementData::calc_num_entities(), moab::ScdElementData::calc_num_entities(), moab::Skinner::classify_2d_boundary(), moab::DualTool::construct_dual_cells(), moab::DualTool::construct_dual_edges(), moab::DualTool::construct_dual_faces(), moab::DualTool::construct_dual_vertices(), moab::HigherOrderFactory::convert_sequence(), moab::SpectralMeshTool::convert_to_coarse(), moab::HigherOrderFactory::copy_mid_face_nodes(), moab::SequenceManager::create_scd_sequence(), moab::ReadGmsh::create_sets(), moab::Skinner::create_side(), moab::ReadNCDF::create_sideset_element(), moab::SequenceManager::create_sweep_sequence(), moab::MeshSet::DIM_FROM_HANDLE(), moab::Core::dimension_from_handle(), moab::Skinner::face_reversed(), moab::ParallelComm::find_existing_entity(), moab::Skinner::find_skin_noadj(), moab::Skinner::find_skin_vertices(), moab::WriteNCDF::gather_mesh_information(), moab::WriteSLAC::gather_mesh_information(), moab::WriteTemplate::gather_mesh_information(), moab::AEntityFactory::get_adjacencies(), moab::get_adjacencies_intersection(), moab::get_adjacencies_union(), moab::MeshTopoUtil::get_bridge_adjacencies(), moab::AEntityFactory::get_element(), moab::AEntityFactory::get_elements(), moab::WriteCCMIO::get_neuset_elems(), moab::WriteSLAC::get_neuset_elems(), moab::WriteTemplate::get_neuset_elems(), moab::VectorSetIterator::get_next_arr(), moab::RangeSetIterator::get_next_by_dimension(), moab::ReadUtil::get_ordered_vertices(), moab::ScdElementData::get_params_connectivity(), moab::SweptElementData::get_params_connectivity(), moab::WriteNCDF::get_sideset_elems(), moab::AEntityFactory::get_up_adjacency_elements(), moab::WriteSLAC::get_valid_sides(), moab::WriteTemplate::get_valid_sides(), moab::WriteNCDF::get_valid_sides(), moab::Core::high_order_node(), HONodeParent(), moab::WriteHDF5::initialize_mesh(), moab::Core::list_entity(), moab::ReadHDF5::load_file_partial(), moab::AEntityFactory::merge_adjust_adjacencies(), moab::Core::merge_entities(), moab::Range::num_of_dimension(), orient_faces_outward(), moab::Tqdcfr::read_block(), moab::Tqdcfr::read_elements(), moab::Tqdcfr::GeomHeader::read_info_header(), moab::AEntityFactory::remove_all_adjacencies(), moab::HigherOrderFactory::remove_mid_face_nodes(), moab::Core::side_number(), SideNumber(), moab::MeshTopoUtil::split_entities_manifold(), moab::ExoIIUtil::static_get_element_type(), SubEntityType(), moab::ReadNCDF::update(), and moab::HigherOrderFactory::zero_mid_face_nodes().

◆ EntityTypeFromName()

EntityType moab::CN::EntityTypeFromName ( const char *  name)
static

given a name, find the corresponding entity type

return a type for the given name

Definition at line 56 of file CN.cpp.

57 {
58  for( EntityType i = MBVERTEX; i < MBMAXTYPE; i++ )
59  {
60  if( 0 == strcmp( name, entityTypeNames[i] ) ) return i;
61  }
62 
63  return MBMAXTYPE;
64 }

References entityTypeNames, MBMAXTYPE, and MBVERTEX.

Referenced by moab::ReadHDF5::load_file_impl(), moab::ReadHDF5::load_file_partial(), moab::ReadHDF5::read_elems(), moab::ReadHDF5::read_node_adj_elems(), and moab::ReadHDF5::read_poly().

◆ EntityTypeName()

const char * moab::CN::EntityTypeName ( const EntityType  this_type)
static

◆ GetBasis()

short int moab::CN::GetBasis ( )
inlinestatic

get the basis of the numbering system

Definition at line 532 of file CN.hpp.

533 {
534  return numberBasis;
535 }

References numberBasis.

◆ getDimPair()

DimensionPair moab::CN::getDimPair ( int  entity_type)
static

Get the dimension pair corresponding to a dimension.

Definition at line 50 of file CN.cpp.

51 {
52  return TypeDimensionMap[entity_type];
53 }

References TypeDimensionMap.

◆ HasMidEdgeNodes()

bool moab::CN::HasMidEdgeNodes ( const EntityType  this_type,
const int  num_verts 
)
inlinestatic

true if entities of a given type and number of nodes indicates mid edge nodes are present.

Parameters
this_typeType of entity for which sub-entity connectivity is being queried
num_vertsNumber of nodes defining entity
Returns
bool Returns true if this_type combined with num_nodes indicates mid-edge nodes are likely

Definition at line 549 of file CN.hpp.

550 {
551  const int bits = HasMidNodes( this_type, num_nodes );
552  return static_cast< bool >( ( bits & ( 1 << 1 ) ) >> 1 );
553 }

References HasMidNodes().

Referenced by moab::HigherOrderFactory::center_node_exist(), moab::Skinner::find_skin_vertices_2D(), and moab::ElementSequence::has_mid_edge_nodes().

◆ HasMidFaceNodes()

bool moab::CN::HasMidFaceNodes ( const EntityType  this_type,
const int  num_verts 
)
inlinestatic

true if entities of a given type and number of nodes indicates mid face nodes are present.

Parameters
this_typeType of entity for which sub-entity connectivity is being queried
num_vertsNumber of nodes defining entity
Returns
bool Returns true if this_type combined with num_nodes indicates mid-face nodes are likely

Definition at line 555 of file CN.hpp.

556 {
557  const int bits = HasMidNodes( this_type, num_nodes );
558  return static_cast< bool >( ( bits & ( 1 << 2 ) ) >> 2 );
559 }

References HasMidNodes().

Referenced by moab::HigherOrderFactory::center_node_exist(), and moab::ElementSequence::has_mid_face_nodes().

◆ HasMidNodes() [1/2]

int moab::CN::HasMidNodes ( const EntityType  this_type,
const int  num_verts 
)
inlinestatic

Same as above, except returns a single integer with the bits, from least significant to most significant set to one if the corresponding mid nodes on sub entities of the least dimension (0) to the highest dimension (3) are present in the elment type.

Definition at line 567 of file CN.hpp.

568 {
569  assert( (unsigned)num_nodes <= (unsigned)MAX_NODES_PER_ELEMENT );
570  return midNodesPerType[this_type][num_nodes];
571 }

References MAX_NODES_PER_ELEMENT, and midNodesPerType.

◆ HasMidNodes() [2/2]

void moab::CN::HasMidNodes ( const EntityType  this_type,
const int  num_verts,
int  mid_nodes[4] 
)
inlinestatic

true if entities of a given type and number of nodes indicates mid edge/face/region nodes are present.

Parameters
this_typeType of entity for which sub-entity connectivity is being queried
num_vertsNumber of nodes defining entity
mid_nodesIf mid_nodes[i], i=1..2 is non-zero, indicates that mid-edge (i=1), mid-face (i=2), and/or mid-region (i=3) nodes are likely

Definition at line 573 of file CN.hpp.

574 {
575  const int bits = HasMidNodes( this_type, num_nodes );
576  mid_nodes[0] = 0;
577  mid_nodes[1] = ( bits & ( 1 << 1 ) ) >> 1;
578  mid_nodes[2] = ( bits & ( 1 << 2 ) ) >> 2;
579  mid_nodes[3] = ( bits & ( 1 << 3 ) ) >> 3;
580 }

Referenced by moab::Skinner::find_skin_vertices_3D(), moab::AEntityFactory::get_down_adjacency_elements(), HasMidEdgeNodes(), HasMidFaceNodes(), HasMidRegionNodes(), moab::Core::high_order_node(), HONodeIndex(), HONodeParent(), moab::ReadNCDF::read_elements(), moab::Tqdcfr::BlockHeader::read_info_header(), SubEntityNodeIndices(), and moab::WriteCCMIO::write_cells_and_faces().

◆ HasMidRegionNodes()

bool moab::CN::HasMidRegionNodes ( const EntityType  this_type,
const int  num_verts 
)
inlinestatic

true if entities of a given type and number of nodes indicates mid region nodes are present.

Parameters
this_typeType of entity for which sub-entity connectivity is being queried
num_vertsNumber of nodes defining entity
Returns
bool Returns true if this_type combined with num_nodes indicates mid-region nodes are likely

Definition at line 561 of file CN.hpp.

562 {
563  const int bits = HasMidNodes( this_type, num_nodes );
564  return static_cast< bool >( ( bits & ( 1 << 3 ) ) >> 3 );
565 }

References HasMidNodes().

Referenced by moab::ElementSequence::has_mid_volume_nodes().

◆ HONodeIndex()

short int moab::CN::HONodeIndex ( const EntityType  this_type,
const int  num_verts,
const int  subfacet_dim,
const int  subfacet_index 
)
static

for an entity of this type with num_verts vertices, and a specified subfacet (dimension and index), return the index of the higher order node for that entity in this entity's connectivity array

for an entity of this type and a specified subfacet (dimension and index), return the index of the higher order node for that entity in this entity's connectivity array

Parameters
this_typeType of entity being queried
num_vertsNumber of vertices for the entity being queried
subfacet_dimDimension of sub-entity being queried
subfacet_indexIndex of sub-entity being queried
Returns
index Index of sub-entity's higher-order node

Definition at line 588 of file CN.cpp.

592 {
593  int i;
594  int has_mids[4];
595  HasMidNodes( this_type, num_verts, has_mids );
596 
597  // if we have no mid nodes on the subfacet_dim, we have no index
598  if( subfacet_index != -1 && !has_mids[subfacet_dim] ) return -1;
599 
600  // put start index at last index (one less than the number of vertices
601  // plus the index basis)
602  int index = VerticesPerEntity( this_type ) - 1 + numberBasis;
603 
604  // for each subfacet dimension less than the target subfacet dim which has mid nodes,
605  // add the number of subfacets of that dimension to the index
606  for( i = 1; i < subfacet_dim; i++ )
607  if( has_mids[i] ) index += NumSubEntities( this_type, i );
608 
609  // now add the index of this subfacet, or one if we're asking about the entity as a whole
610  if( subfacet_index == -1 && has_mids[subfacet_dim] )
611  // want the index of the last ho node on this subfacet
612  index += NumSubEntities( this_type, subfacet_dim );
613 
614  else if( subfacet_index != -1 && has_mids[subfacet_dim] )
615  index += subfacet_index + 1 - numberBasis;
616 
617  // that's it
618  return index;
619 }

References HasMidNodes(), numberBasis, NumSubEntities(), and VerticesPerEntity().

Referenced by moab::Skinner::find_skin_vertices_2D(), moab::AEntityFactory::get_down_adjacency_elements(), and SubEntityNodeIndices().

◆ HONodeParent()

void moab::CN::HONodeParent ( EntityType  elem_type,
int  num_nodes,
int  ho_node_index,
int &  parent_dim,
int &  parent_index 
)
static

given data about an element and a vertex in that element, return the dimension and index of the sub-entity that the vertex resolves. If it does not resolve a sub-entity, either because it's a corner node or it's not in the element, -1 is returned in both return values.

given data about an element and a vertex in that element, return the dimension and index of the sub-entity that the vertex resolves. If it does not resolve a sub-entity, either because it's a corner node or it's not in the element, -1 is returned in both return values

Parameters
elem_typeType of entity being queried
num_nodesThe number of nodes in the element connectivity
ho_node_indexThe position of the HO node in the connectivity list (zero based)
parent_dimDimension of sub-entity high-order node resolves (returned)
parent_indexIndex of sub-entity high-order node resolves (returned)

Definition at line 625 of file CN.cpp.

626 {
627  // begin with error values
628  parent_dim = parent_index = -1;
629 
630  // given the number of verts and the element type, get the hasmidnodes solution
631  int has_mids[4];
632  HasMidNodes( elem_type, num_verts, has_mids );
633 
634  int index = VerticesPerEntity( elem_type ) - 1;
635  const int dim = Dimension( elem_type );
636 
637  // keep a running sum of the ho node indices for this type of element, and stop
638  // when you get to the dimension which has the ho node
639  for( int i = 1; i < dim; i++ )
640  {
641  if( has_mids[i] )
642  {
643  if( ho_index <= index + NumSubEntities( elem_type, i ) )
644  {
645  // the ho_index resolves an entity of dimension i, so set the return values
646  // and break out of the loop
647  parent_dim = i;
648  parent_index = ho_index - index - 1;
649  return;
650  }
651  else
652  {
653  index += NumSubEntities( elem_type, i );
654  }
655  }
656  }
657 
658  // mid region node case
659  if( has_mids[dim] && ho_index == index + 1 )
660  {
661  parent_dim = dim;
662  parent_index = 0;
663  }
664 }

References dim, Dimension(), HasMidNodes(), NumSubEntities(), and VerticesPerEntity().

Referenced by moab::HigherOrderFactory::tag_for_deletion().

◆ NumSubEntities()

◆ OppositeSide()

short int moab::CN::OppositeSide ( const EntityType  parent_type,
const int  child_index,
const int  child_dim,
int &  opposite_index,
int &  opposite_dim 
)
static

return the dimension and index of the opposite side, given parent entity type and child dimension and index. This function is only defined for certain types of parent/child types: (Parent, Child dim->Opposite dim): (Tri, 1->0), (Tri, 0->1), (Quad, 1->1), (Quad, 0->0), (Tet, 2->0), (Tet, 1->1), (Tet, 0->2), (Hex, 2->2), (Hex, 1->1)(diagonally across element), (Hex, 0->0) (diagonally across element) All other parent types and child dimensions return an error.

Parameters
parent_typeThe type of parent element
child_typeThe type of child element
child_indexThe index of the child element
opposite_indexThe index of the opposite element
Returns
status Returns 0 if successful, -1 if not

Definition at line 379 of file CN.cpp.

384 {
385  switch( parent_type )
386  {
387  case MBEDGE:
388  if( 0 != child_dim )
389  return -1;
390  else
391  opposite_index = 1 - child_index;
392  opposite_dim = 0;
393  break;
394 
395  case MBTRI:
396  switch( child_dim )
397  {
398  case 0:
399  opposite_dim = 1;
400  opposite_index = ( child_index + 1 ) % 3;
401  break;
402  case 1:
403  opposite_dim = 0;
404  opposite_index = ( child_index + 2 ) % 3;
405  break;
406  default:
407  return -1;
408  }
409  break;
410 
411  case MBQUAD:
412  switch( child_dim )
413  {
414  case 0:
415  case 1:
416  opposite_dim = child_dim;
417  opposite_index = ( child_index + 2 ) % 4;
418  break;
419  default:
420  return -1;
421  }
422  break;
423 
424  case MBTET:
425  switch( child_dim )
426  {
427  case 0:
428  opposite_dim = 2;
429  opposite_index = ( child_index + 1 ) % 3 + 2 * ( child_index / 3 );
430  break;
431  case 1:
432  opposite_dim = 1;
433  opposite_index = child_index < 3 ? 3 + ( child_index + 2 ) % 3 : ( child_index + 1 ) % 3;
434  break;
435  case 2:
436  opposite_dim = 0;
437  opposite_index = ( child_index + 2 ) % 3 + child_index / 3;
438  break;
439  default:
440  return -1;
441  }
442  break;
443  case MBHEX:
444  opposite_dim = child_dim;
445  switch( child_dim )
446  {
447  case 0:
448  opposite_index = child_index < 4 ? 4 + ( child_index + 2 ) % 4 : ( child_index - 2 ) % 4;
449  break;
450  case 1:
451  opposite_index = 4 * ( 2 - child_index / 4 ) + ( child_index + 2 ) % 4;
452  break;
453  case 2:
454  opposite_index = child_index < 4 ? ( child_index + 2 ) % 4 : 9 - child_index;
455  break;
456  default:
457  return -1;
458  }
459  break;
460 
461  default:
462  return -1;
463  }
464 
465  return 0;
466 }

References MBEDGE, MBHEX, MBQUAD, MBTET, and MBTRI.

Referenced by moab::MeshTopoUtil::opposite_entity().

◆ permuteThis() [1/4]

int moab::CN::permuteThis ( const EntityType  t,
const int  dim,
int *  pvec,
const int  indices_per_ent,
const int  num_entries 
)
inlinestatic

Permute a handle array according to permutation vector set with setPermute; permutation is done in-place.

Permute this vector.

Parameters
tEntityType of handles in pvec
dimDimension of handles in pvec
pvecHandle array being permuted
indices_per_entNumber of indices per entity
num_entriesNumber of entities in pvec

Definition at line 1115 of file CN.cpp.

1116 {
1117  return permute_this( t, dim, pvec, num_indices, num_entries );
1118 }

References dim, moab::permute_this(), and t.

◆ permuteThis() [2/4]

int moab::CN::permuteThis ( const EntityType  t,
const int  dim,
long *  pvec,
const int  indices_per_ent,
const int  num_entries 
)
inlinestatic

Definition at line 1127 of file CN.cpp.

1132 {
1133  return permute_this( t, dim, pvec, num_indices, num_entries );
1134 }

References dim, moab::permute_this(), and t.

◆ permuteThis() [3/4]

int moab::CN::permuteThis ( const EntityType  t,
const int  dim,
unsigned int *  pvec,
const int  indices_per_ent,
const int  num_entries 
)
inlinestatic

Definition at line 1119 of file CN.cpp.

1124 {
1125  return permute_this( t, dim, pvec, num_indices, num_entries );
1126 }

References dim, moab::permute_this(), and t.

◆ permuteThis() [4/4]

int moab::CN::permuteThis ( const EntityType  t,
const int  dim,
void **  pvec,
const int  indices_per_ent,
const int  num_entries 
)
inlinestatic

Definition at line 1135 of file CN.cpp.

1140 {
1141  return permute_this( t, dim, pvec, num_indices, num_entries );
1142 }

References dim, moab::permute_this(), and t.

◆ resetPermutation()

void moab::CN::resetPermutation ( const EntityType  t,
const int  dim 
)
inlinestatic

Reset permutation or reverse permutation vector.

Parameters
tEntityType whose permutation vector is being reset
dimDimension of facets being reset; if -1 is input, all dimensions are reset

Definition at line 606 of file CN.hpp.

607 {
608  if( -1 == dim )
609  {
610  for( unsigned int i = 0; i < 3; i++ )
611  resetPermutation( t, i );
612  return;
613  }
614 
615  for( short unsigned int i = 0; i < MAX_SUB_ENTITIES; i++ )
616  {
617  revPermuteVec[t][dim][i] = permuteVec[t][dim][i] = i;
618  }
619 
621 }

References dim, moab::MAX_SUB_ENTITIES, permuteVec, revPermuteVec, and t.

◆ revPermuteThis() [1/4]

int moab::CN::revPermuteThis ( const EntityType  t,
const int  dim,
int *  pvec,
const int  indices_per_ent,
const int  num_entries 
)
inlinestatic

Reverse permute a handle array according to reverse permutation vector set with setPermute; reverse permutation is done in-place.

Reverse permute this vector.

Parameters
tEntityType of handles in pvec
dimDimension of handles in pvec
pvecHandle array being reverse permuted
indices_per_entNumber of indices per entity
num_entriesNumber of entities in pvec

Definition at line 1145 of file CN.cpp.

1150 {
1151  return rev_permute_this( t, dim, pvec, num_indices, num_entries );
1152 }

References dim, moab::rev_permute_this(), and t.

◆ revPermuteThis() [2/4]

int moab::CN::revPermuteThis ( const EntityType  t,
const int  dim,
long *  pvec,
const int  indices_per_ent,
const int  num_entries 
)
inlinestatic

Definition at line 1161 of file CN.cpp.

1166 {
1167  return rev_permute_this( t, dim, pvec, num_indices, num_entries );
1168 }

References dim, moab::rev_permute_this(), and t.

◆ revPermuteThis() [3/4]

int moab::CN::revPermuteThis ( const EntityType  t,
const int  dim,
unsigned int *  pvec,
const int  indices_per_ent,
const int  num_entries 
)
inlinestatic

Definition at line 1153 of file CN.cpp.

1158 {
1159  return rev_permute_this( t, dim, pvec, num_indices, num_entries );
1160 }

References dim, moab::rev_permute_this(), and t.

◆ revPermuteThis() [4/4]

int moab::CN::revPermuteThis ( const EntityType  t,
const int  dim,
void **  pvec,
const int  indices_per_ent,
const int  num_entries 
)
inlinestatic

Definition at line 1169 of file CN.cpp.

1174 {
1175  return rev_permute_this( t, dim, pvec, num_indices, num_entries );
1176 }

References dim, moab::rev_permute_this(), and t.

◆ SetBasis()

void moab::CN::SetBasis ( const int  in_basis)
static

set the basis of the numbering system

set the basis of the numbering system; may or may not do things besides setting the member variable

Definition at line 44 of file CN.cpp.

45 {
46  numberBasis = in_basis;
47 }

References numberBasis.

◆ setPermutation()

void moab::CN::setPermutation ( const EntityType  t,
const int  dim,
short int *  pvec,
const int  num_entries,
const bool  is_reverse = false 
)
inlinestatic

Set permutation or reverse permutation vector Forward permutation is from CN's numbering into application's ordering; that is, if i is CN's index, pvec[i] is application's index. This function stores the permutation vector for this type and facet dimension, which then is used in calls to permuteThis or revPermuteThis.

Set permutation or reverse permutation vector.

Parameters
tEntityType for which to set permutation
dimDimension of facets whose permutation array is being set
pvecPermutation array
num_entriesNumber of indicies in permutation array
is_reverseArray is reverse permutation

Definition at line 583 of file CN.hpp.

588 {
589  short int *this_vec = permuteVec[t][dim], *that_vec = revPermuteVec[t][dim];
590  if( is_reverse )
591  {
592  this_vec = revPermuteVec[t][dim];
593  that_vec = permuteVec[t][dim];
594  }
595 
596  for( short int i = 0; i < num_entries; i++ )
597  {
598  this_vec[i] = pvec[i];
599  that_vec[pvec[i]] = i;
600  }
601 
602  this_vec[MAX_SUB_ENTITIES] = that_vec[MAX_SUB_ENTITIES] = (short)num_entries;
603 }

References dim, moab::MAX_SUB_ENTITIES, permuteVec, revPermuteVec, and t.

◆ SideNumber() [1/7]

short int moab::CN::SideNumber ( const EntityType  parent_type,
const int *  child_conn_indices,
const int  child_num_verts,
const int  child_dim,
int &  side_number,
int &  sense,
int &  offset 
)
static

return the side index represented in the input sub-entity connectivity

Parameters
parent_typeEntity type of parent entity
child_conn_indicesChild connectivity to query, specified as indices into the connectivity list of the parent.
child_num_vertsNumber of values in child_conn_indices
child_dimDimension of child entity being queried
side_numberSide number of child entity (returned)
senseSense of child entity with respect to order in child_conn (returned)
offsetOffset of child_conn with respect to canonical ordering data (returned)
Returns
status Returns zero if successful, -1 if not

Definition at line 311 of file CN.cpp.

318 {
319  int parent_dim = Dimension( parent_type );
320  int parent_num_verts = VerticesPerEntity( parent_type );
321 
322  // degenerate case (vertex), output == input
323  if( child_dim == 0 )
324  {
325  if( child_num_verts != 1 ) return -1;
326  side_no = *child_conn_indices;
327  sense = offset = 0;
328  }
329 
330  // given a parent and child element, find the corresponding side number
331 
332  // dim_diff should be -1, 0 or 1 (same dimension, one less dimension, two less, resp.)
333  if( child_dim > parent_dim || child_dim < 0 ) return -1;
334 
335  // different types of same dimension won't be the same
336  if( parent_dim == child_dim && parent_num_verts != child_num_verts )
337  {
338  side_no = -1;
339  sense = 0;
340  return 0;
341  }
342 
343  // loop over the sub-elements, comparing to child connectivity
344  int sub_conn_indices[10];
345  for( int i = 0; i < NumSubEntities( parent_type, child_dim ); i++ )
346  {
347  int sub_size = VerticesPerEntity( SubEntityType( parent_type, child_dim, i ) );
348  if( sub_size != child_num_verts ) continue;
349 
350  SubEntityVertexIndices( parent_type, child_dim, i, sub_conn_indices );
351  bool they_match = ConnectivityMatch( child_conn_indices, sub_conn_indices, sub_size, sense, offset );
352  if( they_match )
353  {
354  side_no = i;
355  return 0;
356  }
357  }
358 
359  // if we've gotten here, we don't match
360  side_no = -1;
361 
362  // return value is no success
363  return 1;
364 }

References ConnectivityMatch(), Dimension(), NumSubEntities(), SubEntityType(), SubEntityVertexIndices(), and VerticesPerEntity().

◆ SideNumber() [2/7]

short int moab::CN::SideNumber ( const EntityType  parent_type,
const int *  parent_conn,
const int *  child_conn,
const int  child_num_verts,
const int  child_dim,
int &  side_number,
int &  sense,
int &  offset 
)
static

return the side index represented in the input sub-entity connectivity in the input parent entity connectivity array.

Parameters
parent_connConnectivity of parent entity being queried
parent_typeEntity type of parent entity
child_connConnectivity of child whose index is being queried
child_num_vertsNumber of vertices in child_conn
child_dimDimension of child entity being queried
side_numberSide number of child entity (returned)
senseSense of child entity with respect to order in child_conn (returned)
offsetOffset of child_conn with respect to canonical ordering data (returned)
Returns
status Returns zero if successful, -1 if not

Definition at line 240 of file CN.cpp.

248 {
249  return side_number( parent_conn, parent_type, child_conn, child_num_verts, child_dim, side_no, sense, offset );
250 }

References moab::side_number().

Referenced by moab::Skinner::create_side(), moab::Skinner::face_reversed(), moab::Skinner::find_skin_vertices_2D(), moab::Skinner::find_skin_vertices_3D(), moab::Core::high_order_node(), orient_faces_outward(), moab::GeomTopoTool::restore_topology_from_adjacency(), SphereDecomp::retrieve_subdiv_verts(), moab::Core::side_number(), moab::side_number(), SubEntityNodeIndices(), and moab::WriteCCMIO::write_cells_and_faces().

◆ SideNumber() [3/7]

short int moab::CN::SideNumber ( const EntityType  parent_type,
const long *  parent_conn,
const long *  child_conn,
const int  child_num_verts,
const int  child_dim,
int &  side_number,
int &  sense,
int &  offset 
)
static

Definition at line 263 of file CN.cpp.

271 {
272  return side_number( parent_conn, parent_type, child_conn, child_num_verts, child_dim, side_no, sense, offset );
273 }

References moab::side_number().

◆ SideNumber() [4/7]

short int moab::CN::SideNumber ( const EntityType  parent_type,
const unsigned int *  parent_conn,
const unsigned int *  child_conn,
const int  child_num_verts,
const int  child_dim,
int &  side_number,
int &  sense,
int &  offset 
)
static

Definition at line 252 of file CN.cpp.

260 {
261  return side_number( parent_conn, parent_type, child_conn, child_num_verts, child_dim, side_no, sense, offset );
262 }

References moab::side_number().

◆ SideNumber() [5/7]

short int moab::CN::SideNumber ( const EntityType  parent_type,
const unsigned long *  parent_conn,
const unsigned long *  child_conn,
const int  child_num_verts,
const int  child_dim,
int &  side_number,
int &  sense,
int &  offset 
)
static

Definition at line 274 of file CN.cpp.

282 {
283  return side_number( parent_conn, parent_type, child_conn, child_num_verts, child_dim, side_no, sense, offset );
284 }

References moab::side_number().

◆ SideNumber() [6/7]

short int moab::CN::SideNumber ( const EntityType  parent_type,
const unsigned long long *  parent_conn,
const unsigned long long *  child_conn,
const int  child_num_verts,
const int  child_dim,
int &  side_number,
int &  sense,
int &  offset 
)
static

Definition at line 286 of file CN.cpp.

295 {
296  return side_number( parent_conn, parent_type, child_conn, child_num_verts, child_dim, side_no, sense, offset );
297 }

References moab::side_number().

◆ SideNumber() [7/7]

short int moab::CN::SideNumber ( const EntityType  parent_type,
void *const *  parent_conn,
void *const *  child_conn,
const int  child_num_verts,
const int  child_dim,
int &  side_number,
int &  sense,
int &  offset 
)
static

Definition at line 299 of file CN.cpp.

307 {
308  return side_number( parent_conn, parent_type, child_conn, child_num_verts, child_dim, side_no, sense, offset );
309 }

References moab::side_number().

◆ SubEntityConn()

void moab::CN::SubEntityConn ( const void *  parent_conn,
const EntityType  parent_type,
const int  sub_dimension,
const int  sub_index,
void *  sub_entity_conn,
int &  num_sub_vertices 
)
static

return the vertices of the specified sub entity

Parameters
parent_connConnectivity of parent entity
parent_typeEntity type of parent entity
sub_dimensionDimension of sub-entity being queried
sub_indexIndex of sub-entity being queried
sub_entity_connConnectivity of sub-entity, based on parent_conn and canonical ordering for parent_type
num_sub_verticesNumber of vertices in sub-entity

Definition at line 119 of file CN.cpp.

125 {
126  static int sub_indices[MAX_SUB_ENTITY_VERTICES];
127 
128  SubEntityVertexIndices( parent_type, sub_dimension, sub_index, sub_indices );
129 
130  num_sub_vertices = VerticesPerEntity( SubEntityType( parent_type, sub_dimension, sub_index ) );
131  void** parent_conn_ptr = static_cast< void** >( const_cast< void* >( parent_conn ) );
132  void** sub_conn_ptr = static_cast< void** >( sub_entity_conn );
133  for( int i = 0; i < num_sub_vertices; i++ )
134  sub_conn_ptr[i] = parent_conn_ptr[sub_indices[i]];
135 }

References moab::MAX_SUB_ENTITY_VERTICES, SubEntityType(), SubEntityVertexIndices(), and VerticesPerEntity().

Referenced by moab::WriteCCMIO::write_cells_and_faces().

◆ SubEntityNodeIndices()

void moab::CN::SubEntityNodeIndices ( const EntityType  this_topo,
const int  num_nodes,
const int  sub_dimension,
const int  sub_index,
EntityType &  sub_entity_topo,
int &  num_sub_entity_nodes,
int  sub_entity_conn[] 
)
static

return the node indices of the specified sub-entity.

Parameters
this_topoThe topology of the queried element type
num_nodesThe number of nodes in the queried element type.
sub_dimensionDimension of sub-entity
sub_indexIndex of sub-entity
sub_entity_topo(Output) Topology of requested sub-entity.
num_sub_entity_nodes(Output) Number of nodes in the requested sub-entity.
sub_entity_conn(Output) Connectivity of sub-entity

Definition at line 66 of file CN.cpp.

73 {
74  // If asked for a node, the special case...
75  if( sub_dimension == 0 )
76  {
77  assert( sub_index < num_nodes );
78  subentity_topo = MBVERTEX;
79  num_sub_entity_nodes = 1;
80  sub_entity_conn[0] = sub_index;
81  return;
82  }
83 
84  const int ho_bits = HasMidNodes( this_topo, num_nodes );
85  subentity_topo = SubEntityType( this_topo, sub_dimension, sub_index );
86  num_sub_entity_nodes = VerticesPerEntity( subentity_topo );
87  const short* corners = mConnectivityMap[this_topo][sub_dimension - 1].conn[sub_index];
88  std::copy( corners, corners + num_sub_entity_nodes, sub_entity_conn );
89 
90  int sub_sub_corners[MAX_SUB_ENTITY_VERTICES];
91  int side, sense, offset;
92  for( int dim = 1; dim <= sub_dimension; ++dim )
93  {
94  if( !( ho_bits & ( 1 << dim ) ) ) continue;
95 
96  const short num_mid = NumSubEntities( subentity_topo, dim );
97  for( int i = 0; i < num_mid; ++i )
98  {
99  const EntityType sub_sub_topo = SubEntityType( subentity_topo, dim, i );
100  const int sub_sub_num_vert = VerticesPerEntity( sub_sub_topo );
101  SubEntityVertexIndices( subentity_topo, dim, i, sub_sub_corners );
102 
103  for( int j = 0; j < sub_sub_num_vert; ++j )
104  sub_sub_corners[j] = corners[sub_sub_corners[j]];
105  SideNumber( this_topo, sub_sub_corners, sub_sub_num_vert, dim, side, sense, offset );
106  sub_entity_conn[num_sub_entity_nodes++] = HONodeIndex( this_topo, num_nodes, dim, side );
107  }
108  }
109 }

References moab::CN::ConnMap::conn, dim, HasMidNodes(), HONodeIndex(), moab::MAX_SUB_ENTITY_VERTICES, MBVERTEX, mConnectivityMap, NumSubEntities(), SideNumber(), SubEntityType(), SubEntityVertexIndices(), and VerticesPerEntity().

Referenced by moab::Skinner::classify_2d_boundary(), moab::Skinner::create_side(), moab::ReadNCDF::create_ss_elements(), moab::Skinner::find_skin_noadj(), and moab::Skinner::find_skin_vertices_3D().

◆ SubEntityType()

EntityType moab::CN::SubEntityType ( const EntityType  this_type,
const int  sub_dimension,
const int  index 
)
static

return the type of a particular sub-entity.

Parameters
this_typeType of entity for which sub-entity type is being queried
sub_dimensionTopological dimension of sub-entity whose type is being queried
indexIndex of sub-entity whose type is being queried
Returns
type Entity type of sub-entity with specified dimension and index

Definition at line 1085 of file CN.cpp.

1086 {
1087  return ( !sub_dimension ? MBVERTEX
1088  : ( Dimension( this_type ) == sub_dimension && 0 == index
1089  ? this_type
1090  : mConnectivityMap[this_type][sub_dimension - 1].target_type[index] ) );
1091 }

References Dimension(), MBVERTEX, and mConnectivityMap.

Referenced by moab::GeomUtil::box_linear_elem_overlap(), moab::MeshTopoUtil::get_bridge_adjacencies(), SideNumber(), SubEntityConn(), and SubEntityNodeIndices().

◆ SubEntityVertexIndices() [1/2]

const short * moab::CN::SubEntityVertexIndices ( const EntityType  this_type,
const int  sub_dimension,
const int  sub_index,
EntityType &  sub_type,
int &  num_sub_ent_vertices 
)
static

return the vertex indices of the specified sub-entity.

Parameters
this_typeType of entity for which sub-entity connectivity is being queried
sub_dimensionDimension of sub-entity
sub_indexIndex of sub-entity
num_sub_ent_verticesthe number of vertices in the sub-entity

Definition at line 1093 of file CN.cpp.

1098 {
1099  if( sub_dimension == 0 )
1100  {
1101  n = 1;
1102  sub_type = MBVERTEX;
1103  return increasingInts + index;
1104  }
1105  else
1106  {
1107  const CN::ConnMap& map = mConnectivityMap[this_type][sub_dimension - 1];
1108  sub_type = map.target_type[index];
1109  n = map.num_corners_per_sub_element[index];
1110  return map.conn[index];
1111  }
1112 }

References moab::CN::ConnMap::conn, increasingInts, MBVERTEX, mConnectivityMap, moab::CN::ConnMap::num_corners_per_sub_element, and moab::CN::ConnMap::target_type.

◆ SubEntityVertexIndices() [2/2]

void moab::CN::SubEntityVertexIndices ( const EntityType  this_type,
const int  sub_dimension,
const int  sub_index,
int  sub_entity_conn[] 
)
inlinestatic

return the vertex indices of the specified sub-entity.

return the connectivity of the specified sub-entity.

Parameters
this_typeType of entity for which sub-entity connectivity is being queried
sub_dimensionDimension of sub-entity
sub_indexIndex of sub-entity
sub_entity_connConnectivity of sub-entity (returned to calling function)

Definition at line 538 of file CN.hpp.

542 {
543  EntityType type;
544  int n;
545  const short* indices = SubEntityVertexIndices( this_type, sub_dimension, index, type, n );
546  std::copy( indices, indices + n, sub_entity_conn );
547 }

Referenced by moab::GeomUtil::box_linear_elem_overlap(), moab::Skinner::classify_2d_boundary(), moab::Skinner::find_skin_vertices_3D(), gather_set_stats(), moab::MeshTopoUtil::get_bridge_adjacencies(), moab::ReadUtil::get_ordered_vertices(), min_edge_length(), SideNumber(), SubEntityConn(), and SubEntityNodeIndices().

◆ SwitchBasis()

static void moab::CN::SwitchBasis ( const int  old_basis,
const int  new_basis 
)
staticprivate

switch the basis

◆ VerticesPerEntity()

short int moab::CN::VerticesPerEntity ( const EntityType  t)
static

return the number of (corner) vertices contained in the specified type.

Definition at line 1071 of file CN.cpp.

1072 {
1073  return ( MBVERTEX == t
1074  ? (short int)1
1075  : mConnectivityMap[t][mConnectivityMap[t][0].topo_dimension - 1].num_corners_per_sub_element[0] );
1076 }

References MBVERTEX, mConnectivityMap, and t.

Referenced by moab::Skinner::add_adjacency(), moab::HigherOrderFactory::add_mid_edge_nodes(), moab::HigherOrderFactory::add_mid_face_nodes(), moab::HigherOrderFactory::add_mid_volume_nodes(), moab::GeomUtil::box_linear_elem_overlap(), moab::HigherOrderFactory::center_node_exist(), moab::AEntityFactory::check_equiv_entities(), moab::Skinner::classify_2d_boundary(), SphereDecomp::compute_nodes(), moab::HigherOrderFactory::convert_sequence(), moab::HigherOrderFactory::copy_corner_nodes(), moab::HigherOrderFactory::copy_mid_edge_nodes(), moab::HigherOrderFactory::copy_mid_face_nodes(), moab::HigherOrderFactory::copy_mid_volume_nodes(), moab::Core::create_element(), moab::ReadIDEAS::create_elements(), moab::Skinner::create_side(), moab::ReadNCDF::create_sideset_element(), moab::AEntityFactory::entities_equivalent(), moab::Skinner::face_reversed(), moab::Skinner::find_match(), moab::Skinner::find_skin_noadj(), moab::Skinner::find_skin_vertices_3D(), moab::WriteVtk::gather_mesh(), moab::MeshTopoUtil::get_bridge_adjacencies(), moab::UnstructuredElemSeq::get_connectivity(), moab::Core::get_connectivity_by_type(), moab::ReadUtil::get_ordered_vertices(), moab::WriteNCDF::get_valid_sides(), moab::Core::high_order_node(), HONodeIndex(), HONodeParent(), moab::HigherOrderFactory::initialize_map(), moab::ReadSms::load_file_impl(), moab::WriteGMV::local_write_mesh(), NumSubEntities(), moab::Tqdcfr::read_block(), moab::ReadNASTRAN::read_element(), moab::Tqdcfr::BlockHeader::read_info_header(), moab::HigherOrderFactory::remove_mid_edge_nodes(), moab::HigherOrderFactory::remove_mid_face_nodes(), moab::HigherOrderFactory::remove_mid_volume_nodes(), moab::side_number(), SideNumber(), moab::ExoIIUtil::static_get_element_type(), SubEntityConn(), SubEntityNodeIndices(), moab::WriteVtk::write_elems(), moab::HigherOrderFactory::zero_mid_edge_nodes(), moab::HigherOrderFactory::zero_mid_face_nodes(), and moab::HigherOrderFactory::zero_mid_volume_nodes().

Member Data Documentation

◆ entityTypeNames

const char * moab::CN::entityTypeNames
staticprivate
Initial value:
= { "Vertex", "Edge", "Tri", "Quad", "Polygon", "Tet", "Pyramid",
"Prism", "Knife", "Hex", "Polyhedron", "EntitySet", "MaxType" }

entity names

Definition at line 60 of file CN.hpp.

Referenced by EntityTypeFromName(), and EntityTypeName().

◆ increasingInts

short moab::CN::increasingInts
staticprivate
Initial value:
= { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19,
20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39 }

Definition at line 71 of file CN.hpp.

Referenced by SubEntityVertexIndices().

◆ mConnectivityMap

◆ midNodesPerType

const unsigned char moab::CN::midNodesPerType
static
Initial value:
= {
{ 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, E, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, 0, F, 0, E, E | F, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, 0, 0, F, 0, 0, E, E | F, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, 0, 0, R, 0, 0, F, F | R, E, E | R, 0, 0, E | F, E | F | R, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, 0, 0, 0, R, 0, 0, 0, F, F | R, 0, E, E | R, 0, 0, 0, E | F, E | F | R, 0, 0, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, 0, 0, 0, 0, R, 0, 0, 0, F, F | R, 0, 0, E, E | R, 0, 0, 0, E | F, E | F | R, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, 0, 0, 0, 0, 0, R, 0, 0, 0, F, F | R, 0, 0, 0, E, E | R, 0, 0, 0, E | F, E | F | R, 0, 0, 0, 0 },
{ 0, 0, 0, 0, 0, 0, 0, 0, 0, R, 0, 0, 0, 0, F, F | R, 0, 0, 0, 0, E, E | R, 0, 0, 0, 0, E | F, E | F | R },
{ 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
{ 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
}

Definition at line 142 of file CN.hpp.

Referenced by HasMidNodes().

◆ mUpConnMap

const CN::UpConnMap moab::CN::mUpConnMap
static

Definition at line 139 of file CN.hpp.

◆ numberBasis

short int moab::CN::numberBasis = 0
staticprivate

the basis of the numbering system (normally 0 or 1, 0 by default)

Definition at line 66 of file CN.hpp.

Referenced by GetBasis(), HONodeIndex(), and SetBasis().

◆ permuteVec

short int moab::CN::permuteVec
static

Permutation and reverse permutation vectors.

Definition at line 145 of file CN.hpp.

Referenced by moab::permute_this(), resetPermutation(), and setPermutation().

◆ revPermuteVec

short int moab::CN::revPermuteVec
static

Definition at line 146 of file CN.hpp.

Referenced by resetPermutation(), moab::rev_permute_this(), and setPermutation().

◆ TypeDimensionMap


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