moab
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>

List of all members.

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 short int GetBasis ()
 get the basis of the numbering system
static void SetBasis (const int in_basis)
 set the basis of the numbering system
static const char * EntityTypeName (const EntityType this_type)
 return the string type name for this type
static EntityType EntityTypeFromName (const char *name)
 given a name, find the corresponding entity type
static short int Dimension (const EntityType t)
 return the topological entity dimension
static short int VerticesPerEntity (const EntityType t)
 return the number of (corner) vertices contained in the specified type.
static short int NumSubEntities (const EntityType t, const int d)
 return the number of sub-entities bounding the entity.
static EntityType SubEntityType (const EntityType this_type, const int sub_dimension, const int index)
 return the type of a particular sub-entity.
static void SubEntityVertexIndices (const EntityType this_type, const int sub_dimension, const int sub_index, int sub_entity_conn[])
 return the connectivity of the specified sub-entity.
static const short * SubEntityVertexIndices (const EntityType this_type, const int sub_dimension, const int sub_index, EntityType &sub_type, int &num_sub_ent_vertices)
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[])
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)
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)
 given an entity and a target dimension & side number, get that entity
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)
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, 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)
static short int OppositeSide (const EntityType parent_type, const int child_index, const int child_dim, int &opposite_index, int &opposite_dim)
static bool ConnectivityMatch (const int *conn1, const int *conn2, const int num_vertices, int &direct, int &offset)
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 (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.
static void resetPermutation (const EntityType t, const int dim)
 Reset permutation or reverse permutation vector.
static int permuteThis (const EntityType t, const int dim, int *pvec, const int indices_per_ent, const int num_entries)
 Permute this vector.
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 this vector.
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)
static bool HasMidFaceNodes (const EntityType this_type, const int num_verts)
static bool HasMidRegionNodes (const EntityType this_type, const int num_verts)
static void HasMidNodes (const EntityType this_type, const int num_verts, int mid_nodes[4])
static int HasMidNodes (const EntityType this_type, const int num_verts)
static void HONodeParent (EntityType elem_type, int num_nodes, int ho_node_index, int &parent_dim, int &parent_index)
static short int HONodeIndex (const EntityType this_type, const int num_verts, const int subfacet_dim, const int subfacet_index)

Static Public Attributes

static const ConnMap mConnectivityMap [MBMAXTYPE][3]
static const UpConnMap mUpConnMap [MBMAXTYPE][4][4]
static 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.
static short int revPermuteVec [MBMAXTYPE][3][MAX_SUB_ENTITIES+1]
static const DimensionPair TypeDimensionMap []

Private Member Functions

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

Static Private Member Functions

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

Static Private Attributes

static const char * entityTypeNames []
 entity names
static short int numberBasis = 0
 the basis of the numbering system (normally 0 or 1, 0 by default)
static 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 Coroporation, 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 52 of file CN.hpp.


Member Enumeration Documentation

anonymous enum
Enumerator:
MAX_NODES_PER_ELEMENT 

Definition at line 72 of file CN.hpp.

anonymous enum
Enumerator:
MID_EDGE_BIT 
MID_FACE_BIT 
MID_REGION_BIT 

Definition at line 73 of file CN.hpp.

       { MID_EDGE_BIT   = 1<<1,
         MID_FACE_BIT   = 1<<2,
         MID_REGION_BIT = 1<<3 };
anonymous enum

enum used to specify operation type

Enumerator:
INTERSECT 
UNION 

Definition at line 78 of file CN.hpp.

{INTERSECT = 0, UNION};

Constructor & Destructor Documentation

moab::CN::CN ( ) [private]

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


Member Function Documentation

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]

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

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.

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 144 of file CN.cpp.

{
    // first get all the vertex indices
  std::vector<int> tmp_indices;
  const int* it1 = source_indices;

  assert(source_dim >= 0 && source_dim <= 3 &&
         target_dim >= 0 && target_dim <= 3 &&
           // make sure we're not stepping off the end of the array; 
         ((source_dim > 0 && 
           *it1 < mConnectivityMap[this_type][source_dim-1].num_sub_elements) ||
          (source_dim == 0 && 
           *it1 < mConnectivityMap[this_type][Dimension(this_type)-1].num_corners_per_sub_element[0])) && 
         *it1 >= 0);


#define MUC CN::mUpConnMap[this_type][source_dim][target_dim]

    // if we're looking for the vertices of a single side, return them in
    // the canonical ordering; otherwise, return them in sorted order
  if (num_source_indices == 1 && 0 == target_dim && source_dim != target_dim) {

      // element of mConnectivityMap should be for this type and for one
      // less than source_dim, which should give the connectivity of that sub element
    const ConnMap &cm = mConnectivityMap[this_type][source_dim-1];
    std::copy(cm.conn[source_indices[0]],
              cm.conn[source_indices[0]]+cm.num_corners_per_sub_element[source_indices[0]],
              std::back_inserter(index_list));
    return 0;
  }
              
    // now go through source indices, folding adjacencies into target list
  for (it1 = source_indices; it1 != source_indices+num_source_indices; it1++) {
      // *it1 is the side index
      // at start of iteration, index_list has the target list

      // if a union, or first iteration and index list was empty, copy the list
    if (operation_type == CN::UNION || 
        (it1 == source_indices && index_list.empty())) {
      std::copy(MUC.targets_per_source_element[*it1],
                MUC.targets_per_source_element[*it1]+
                MUC.num_targets_per_source_element[*it1],
                std::back_inserter(index_list));
    }
    else {
        // else we're intersecting, and have a non-empty list; intersect with this target list
      tmp_indices.clear();
      for (int i = MUC.num_targets_per_source_element[*it1]-1; i>= 0; i--)
        if (std::find(index_list.begin(), index_list.end(), MUC.targets_per_source_element[*it1][i]) !=
            index_list.end())
          tmp_indices.push_back(MUC.targets_per_source_element[*it1][i]);
//      std::set_intersection(MUC.targets_per_source_element[*it1],
//                            MUC.targets_per_source_element[*it1]+
//                            MUC.num_targets_per_source_element[*it1],
//                            index_list.begin(), index_list.end(),
//                            std::back_inserter(tmp_indices));
      index_list.swap(tmp_indices);

        // if we're at this point and the list is empty, the intersection will be NULL;
        // return if so
      if (index_list.empty()) return 0;
    }
  }
  
  if (operation_type == CN::UNION && num_source_indices != 1) {
      // need to sort then unique the list
    std::sort(index_list.begin(), index_list.end());
    index_list.erase(std::unique(index_list.begin(), index_list.end()), 
                     index_list.end());
  }
  
  return 0;
}
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 522 of file CN.cpp.

{
  return connectivity_match<int>(conn1_i, conn2_i, num_vertices, direct, offset );
}
bool moab::CN::ConnectivityMatch ( const unsigned int *  conn1,
const unsigned int *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
) [static]

Definition at line 530 of file CN.cpp.

{
  return connectivity_match<unsigned int>(conn1_i, conn2_i, num_vertices, direct, offset );
}
bool moab::CN::ConnectivityMatch ( const long *  conn1,
const long *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
) [static]

Definition at line 538 of file CN.cpp.

{
  return connectivity_match<long>(conn1_i, conn2_i, num_vertices, direct, offset );
}
bool moab::CN::ConnectivityMatch ( const unsigned long *  conn1,
const unsigned long *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
) [static]

Definition at line 546 of file CN.cpp.

{
  return connectivity_match<unsigned long>(conn1_i, conn2_i, num_vertices, direct, offset );
}
bool moab::CN::ConnectivityMatch ( void *const *  conn1,
void *const *  conn2,
const int  num_vertices,
int &  direct,
int &  offset 
) [static]

Definition at line 554 of file CN.cpp.

{
  return connectivity_match<void*>(conn1_i, conn2_i, num_vertices, direct, offset );
}
short int moab::CN::Dimension ( const EntityType  t) [inline, static]

return the topological entity dimension

Definition at line 471 of file CN.hpp.

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 65 of file CN.cpp.

{
  for (EntityType i = MBVERTEX; i < MBMAXTYPE; i++) {
    if (0 == strcmp(name, entityTypeNames[i]))
      return i;
  }
  
  return MBMAXTYPE;
}
const char * moab::CN::EntityTypeName ( const EntityType  this_type) [inline, static]

return the string type name for this type

Examples:
GetEntities.cpp, and StructuredMeshSimple.cpp.

Definition at line 466 of file CN.hpp.

{
  return entityTypeNames[this_type];
}
short int moab::CN::GetBasis ( ) [inline, static]

get the basis of the numbering system

Definition at line 464 of file CN.hpp.

{return numberBasis;}
bool moab::CN::HasMidEdgeNodes ( const EntityType  this_type,
const int  num_verts 
) [inline, static]

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 529 of file CN.hpp.

{
  const int bits = HasMidNodes( this_type, num_nodes );
  return static_cast<bool>( (bits & (1<<1)) >> 1 );
}
bool moab::CN::HasMidFaceNodes ( const EntityType  this_type,
const int  num_verts 
) [inline, static]

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 536 of file CN.hpp.

{
  const int bits = HasMidNodes( this_type, num_nodes );
  return static_cast<bool>( (bits & (1<<2)) >> 2 );
}
void moab::CN::HasMidNodes ( const EntityType  this_type,
const int  num_verts,
int  mid_nodes[4] 
) [inline, static]

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 557 of file CN.hpp.

{
  const int bits = HasMidNodes( this_type, num_nodes );
  mid_nodes[0] = 0;
  mid_nodes[1] = (bits & (1<<1)) >> 1;
  mid_nodes[2] = (bits & (1<<2)) >> 2;
  mid_nodes[3] = (bits & (1<<3)) >> 3;
}
int moab::CN::HasMidNodes ( const EntityType  this_type,
const int  num_verts 
) [inline, static]

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 550 of file CN.hpp.

{
  assert( (unsigned)num_nodes <= (unsigned)MAX_NODES_PER_ELEMENT );
  return midNodesPerType[this_type][num_nodes];
}
bool moab::CN::HasMidRegionNodes ( const EntityType  this_type,
const int  num_verts 
) [inline, static]

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 543 of file CN.hpp.

{
  const int bits = HasMidNodes( this_type, num_nodes );
  return static_cast<bool>( (bits & (1<<3)) >> 3 );
}
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

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

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

Definition at line 566 of file CN.cpp.

{
  int i;
  int has_mids[4];
  HasMidNodes(this_type, num_verts, has_mids);

    // if we have no mid nodes on the subfacet_dim, we have no index
  if (subfacet_index != -1 && !has_mids[subfacet_dim]) return -1;

    // put start index at last index (one less than the number of vertices 
    // plus the index basis)
  int index = VerticesPerEntity(this_type) - 1 + numberBasis;

    // for each subfacet dimension less than the target subfacet dim which has mid nodes, 
    // add the number of subfacets of that dimension to the index
  for (i = 1; i < subfacet_dim; i++)
    if (has_mids[i]) index += NumSubEntities(this_type, i);
    

    // now add the index of this subfacet, or one if we're asking about the entity as a whole
  if (subfacet_index == -1 && has_mids[subfacet_dim])
      // want the index of the last ho node on this subfacet
    index += NumSubEntities(this_type, subfacet_dim);
  
  else if (subfacet_index != -1 && has_mids[subfacet_dim])
    index += subfacet_index + 1 - numberBasis;

    // that's it
  return index;
}
void moab::CN::HONodeParent ( EntityType  elem_type,
int  num_verts,
int  ho_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.

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)

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

Definition at line 602 of file CN.cpp.

{
    // begin with error values
  parent_dim = parent_index = -1;
     
    // given the number of verts and the element type, get the hasmidnodes solution
  int has_mids[4];
  HasMidNodes(elem_type, num_verts, has_mids);

  int index = VerticesPerEntity(elem_type)-1;
  const int dim = Dimension(elem_type);

    // keep a running sum of the ho node indices for this type of element, and stop
    // when you get to the dimension which has the ho node
  for (int i = 1; i < dim; i++) {
    if (has_mids[i]) {
      if (ho_index <= index + NumSubEntities(elem_type, i)) {
          // the ho_index resolves an entity of dimension i, so set the return values
          // and break out of the loop
        parent_dim = i;
        parent_index = ho_index - index - 1;
        return;
      }
      else {
        index += NumSubEntities(elem_type, i);
      } 
    }
  }
  
    // mid region node case  
  if( has_mids[dim] && ho_index == index+1 ) {
    parent_dim = dim;
    parent_index = 0; 
  }
}
short int moab::CN::NumSubEntities ( const EntityType  t,
const int  d 
) [inline, static]

return the number of sub-entities bounding the entity.

Definition at line 481 of file CN.hpp.

{
  return (t != MBVERTEX && d > 0 ? mConnectivityMap[t][d-1].num_sub_elements :
          (d ? (short int) -1 : VerticesPerEntity(t)));
}
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 360 of file CN.cpp.

{
  switch (parent_type) {
    case MBEDGE:
        if (0 != child_dim) return -1;
        else opposite_index = 1-child_index;
        opposite_dim = 0;
        break;
      
    case MBTRI:
        switch (child_dim) {
          case 0:
              opposite_dim = 1;
              opposite_index = (child_index+1)%3;
              break;
          case 1:
              opposite_dim = 0;
              opposite_index = (child_index+2)%3;
              break;
          default:
              return -1;
        }
        break;

    case MBQUAD:
        switch (child_dim) {
          case 0:
          case 1:
              opposite_dim = child_dim;
              opposite_index = (child_index+2)%4;
              break;
          default:
              return -1;
        }
        break;
      
    case MBTET:
        switch (child_dim) {
          case 0:
              opposite_dim = 2;
              opposite_index = (child_index+1)%3 + 2*(child_index/3);
              break;
          case 1:
              opposite_dim = 1;
              opposite_index = child_index < 3 
                             ? 3 + (child_index + 2)%3
                             : (child_index + 1)%3;
              break;
          case 2:
              opposite_dim = 0;
              opposite_index = (child_index+2)%3 + child_index/3;
              break;
          default:
              return -1;
        }
        break;
    case MBHEX:
      opposite_dim = child_dim;
      switch (child_dim) {
        case 0:
          opposite_index = child_index < 4 
                         ? 4 + (child_index + 2) % 4
                         : (child_index - 2) % 4;
          break;
        case 1:
          opposite_index = 4*(2-child_index/4) + (child_index+2)%4;
          break;
        case 2:
          opposite_index = child_index < 4 
                         ? (child_index + 2) % 4
                         : 9 - child_index;
          break;
        default:
          return -1;
      }
      break;
        
      
    default:
        return -1;
  }
  
  return 0;
}
int moab::CN::permuteThis ( const EntityType  t,
const int  dim,
int *  pvec,
const int  indices_per_ent,
const int  num_entries 
) [inline, static]

Permute this vector.

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

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 1025 of file CN.cpp.

{return permute_this(t, dim, pvec, num_indices, num_entries);}
int moab::CN::permuteThis ( const EntityType  t,
const int  dim,
unsigned int *  pvec,
const int  indices_per_ent,
const int  num_entries 
) [inline, static]

Definition at line 1028 of file CN.cpp.

{return permute_this(t, dim, pvec, num_indices, num_entries);}
int moab::CN::permuteThis ( const EntityType  t,
const int  dim,
long *  pvec,
const int  indices_per_ent,
const int  num_entries 
) [inline, static]

Definition at line 1031 of file CN.cpp.

{return permute_this(t, dim, pvec, num_indices, num_entries);}
int moab::CN::permuteThis ( const EntityType  t,
const int  dim,
void **  pvec,
const int  indices_per_ent,
const int  num_entries 
) [inline, static]

Definition at line 1034 of file CN.cpp.

{return permute_this(t, dim, pvec, num_indices, num_entries);}
void moab::CN::resetPermutation ( const EntityType  t,
const int  dim 
) [inline, static]

Reset permutation or reverse permutation vector.

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 586 of file CN.hpp.

{
  if (-1 == dim) {
    for (unsigned int i = 0; i < 3; i++) resetPermutation(t, i);
    return;
  }
  
  for (short unsigned int i = 0; i < MAX_SUB_ENTITIES; i++) {
    revPermuteVec[t][dim][i] = permuteVec[t][dim][i] = i;
  }
  
  revPermuteVec[t][dim][MAX_SUB_ENTITIES] = 
    permuteVec[t][dim][MAX_SUB_ENTITIES] = MAX_SUB_ENTITIES+1;
}
int moab::CN::revPermuteThis ( const EntityType  t,
const int  dim,
int *  pvec,
const int  indices_per_ent,
const int  num_entries 
) [inline, static]

Reverse permute this vector.

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

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 1039 of file CN.cpp.

{return rev_permute_this(t, dim, pvec, num_indices, num_entries);}
int moab::CN::revPermuteThis ( const EntityType  t,
const int  dim,
unsigned int *  pvec,
const int  indices_per_ent,
const int  num_entries 
) [inline, static]

Definition at line 1042 of file CN.cpp.

{return rev_permute_this(t, dim, pvec, num_indices, num_entries);}
int moab::CN::revPermuteThis ( const EntityType  t,
const int  dim,
long *  pvec,
const int  indices_per_ent,
const int  num_entries 
) [inline, static]

Definition at line 1045 of file CN.cpp.

{return rev_permute_this(t, dim, pvec, num_indices, num_entries);}
int moab::CN::revPermuteThis ( const EntityType  t,
const int  dim,
void **  pvec,
const int  indices_per_ent,
const int  num_entries 
) [inline, static]

Definition at line 1048 of file CN.cpp.

{return rev_permute_this(t, dim, pvec, num_indices, num_entries);}
void moab::CN::SetBasis ( const int  in_basis) [static]

set the basis of the numbering system

member variable

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

Definition at line 59 of file CN.cpp.

{
  numberBasis = in_basis;
}
void moab::CN::setPermutation ( const EntityType  t,
const int  dim,
short int *  pvec,
const int  num_entries,
const bool  is_reverse = false 
) [inline, static]

Set permutation or reverse permutation vector.

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.

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 568 of file CN.hpp.

{
  short int *this_vec = permuteVec[t][dim], *that_vec = revPermuteVec[t][dim];
  if (is_reverse) {
    this_vec = revPermuteVec[t][dim];
    that_vec = permuteVec[t][dim];
  }

  for (short int i = 0; i < num_entries; i++) {
    this_vec[i] = pvec[i];
    that_vec[pvec[i]] = i;
  }

  this_vec[MAX_SUB_ENTITIES] = that_vec[MAX_SUB_ENTITIES] = (short)num_entries;
}
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 248 of file CN.cpp.

{
  return side_number(parent_conn, parent_type, child_conn, child_num_verts,
                     child_dim, side_no, sense, offset);
}
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 257 of file CN.cpp.

{
  return side_number(parent_conn, parent_type, child_conn, child_num_verts,
                     child_dim, side_no, sense, offset);
}
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 265 of file CN.cpp.

{
  return side_number(parent_conn, parent_type, child_conn, child_num_verts,
                     child_dim, side_no, sense, offset);
}
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 273 of file CN.cpp.

{
  return side_number(parent_conn, parent_type, child_conn, child_num_verts,
                     child_dim, side_no, sense, offset);
}
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 281 of file CN.cpp.

{
  return side_number(parent_conn, parent_type, child_conn, child_num_verts,
                     child_dim, side_no, sense, offset);
}
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 290 of file CN.cpp.

{
  int parent_dim = Dimension(parent_type);
  int parent_num_verts = VerticesPerEntity(parent_type);

    // degenerate case (vertex), output == input
  if (child_dim == 0) {
    if (child_num_verts != 1)
      return -1;
    side_no = *child_conn_indices;
    sense = offset = 0;
  }
    
    // given a parent and child element, find the corresponding side number

    // dim_diff should be -1, 0 or 1 (same dimension, one less dimension, two less, resp.)
  if (child_dim > parent_dim || child_dim < 0)
    return -1;

    // different types of same dimension won't be the same
  if (parent_dim == child_dim &&
      parent_num_verts != child_num_verts) {
    side_no = -1;
    sense = 0;
    return 0;
  }

    // loop over the sub-elements, comparing to child connectivity
  int sub_conn_indices[10];
  for (int i = 0; i < NumSubEntities(parent_type, child_dim); i++) {
    int sub_size = VerticesPerEntity(SubEntityType(parent_type, child_dim, i));
    if (sub_size != child_num_verts) 
      continue;

    SubEntityVertexIndices(parent_type, child_dim, i, sub_conn_indices);
    bool they_match = ConnectivityMatch(child_conn_indices, 
                                        sub_conn_indices, sub_size, 
                                        sense, offset);
    if (they_match) {
      side_no = i;
      return 0;
    }
  }

    // if we've gotten here, we don't match
  side_no = -1;

    // return value is still success, we didn't have any fatal errors or anything
  return 0;
}
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 127 of file CN.cpp.

{
  static int sub_indices[MAX_SUB_ENTITY_VERTICES];
  
  SubEntityVertexIndices(parent_type, sub_dimension, sub_index, sub_indices);
  
  num_sub_vertices = VerticesPerEntity(SubEntityType(parent_type, sub_dimension, sub_index));
  void **parent_conn_ptr = static_cast<void **>(const_cast<void *>(parent_conn));
  void **sub_conn_ptr = static_cast<void **>(sub_entity_conn);
  for (int i = 0; i < num_sub_vertices; i++)
    sub_conn_ptr[i] = parent_conn_ptr[sub_indices[i]];
}
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 75 of file CN.cpp.

{
    // If asked for a node, the special case...
  if (sub_dimension == 0) {
    assert( sub_index < num_nodes );
    subentity_topo = MBVERTEX;
    num_sub_entity_nodes = 1;
    sub_entity_conn[0] = sub_index;
    return;
  }
  
  const int ho_bits = HasMidNodes( this_topo, num_nodes );
  subentity_topo = SubEntityType(this_topo, sub_dimension, sub_index);
  num_sub_entity_nodes = VerticesPerEntity(subentity_topo);
  const short* corners = mConnectivityMap[this_topo][sub_dimension-1].conn[sub_index];
  std::copy( corners, corners+num_sub_entity_nodes, sub_entity_conn );
  
  int sub_sub_corners[MAX_SUB_ENTITY_VERTICES];
  int side, sense, offset;
  for (int dim = 1; dim <= sub_dimension; ++dim) {
    if (!(ho_bits & (1<<dim)))
      continue;
    
    const short num_mid = NumSubEntities( subentity_topo, dim );
    for (int i = 0; i < num_mid; ++i) {
      const EntityType sub_sub_topo = SubEntityType( subentity_topo, dim, i );
      const int sub_sub_num_vert = VerticesPerEntity( sub_sub_topo );
      SubEntityVertexIndices( subentity_topo, dim, i, sub_sub_corners );

      for (int j = 0; j < sub_sub_num_vert; ++j)
        sub_sub_corners[j] = corners[sub_sub_corners[j]];
      SideNumber( this_topo, sub_sub_corners, sub_sub_num_vert, dim, side, sense, offset );
      sub_entity_conn[num_sub_entity_nodes++] = HONodeIndex( this_topo, num_nodes, dim, side );
    }
  }
}
EntityType moab::CN::SubEntityType ( const EntityType  this_type,
const int  sub_dimension,
const int  index 
) [inline, static]

return the type of a particular sub-entity.

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 488 of file CN.hpp.

{
  
  return (!sub_dimension ? MBVERTEX : 
          (Dimension(this_type) == sub_dimension && 0 == index ? this_type :
          mConnectivityMap[this_type][sub_dimension-1].target_type[index]));
}
void moab::CN::SubEntityVertexIndices ( const EntityType  this_type,
const int  sub_dimension,
const int  sub_index,
int  sub_entity_conn[] 
) [inline, static]

return the connectivity of the specified sub-entity.

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
sub_entity_connConnectivity of sub-entity (returned to calling function)

Definition at line 518 of file CN.hpp.

{
  EntityType type;
  int n;
  const short* indices = SubEntityVertexIndices( this_type, sub_dimension, index, type, n );
  std::copy( indices, indices+n, sub_entity_conn );
}
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 
) [inline, 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 498 of file CN.hpp.

{
  if (sub_dimension == 0) {
    n = 1;
    sub_type = MBVERTEX;
    return increasingInts + index;
  }
  else {
    const CN::ConnMap& map = mConnectivityMap[this_type][sub_dimension-1];
    sub_type = map.target_type[index];
    n = map.num_corners_per_sub_element[index];
    return map.conn[index];
  }
}
static void moab::CN::SwitchBasis ( const int  old_basis,
const int  new_basis 
) [static, private]

switch the basis

short int moab::CN::VerticesPerEntity ( const EntityType  t) [inline, static]

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

Definition at line 476 of file CN.hpp.

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

Member Data Documentation

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

entity names

Definition at line 57 of file CN.hpp.

short moab::CN::increasingInts [static, private]
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 68 of file CN.hpp.

Definition at line 109 of file CN.hpp.

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 129 of file CN.hpp.

Definition at line 126 of file CN.hpp.

short int moab::CN::numberBasis = 0 [static, private]

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

Definition at line 63 of file CN.hpp.

short int moab::CN::permuteVec[MBMAXTYPE][3][MAX_SUB_ENTITIES+1] [static]

Permutation and reverse permutation vectors.

Definition at line 132 of file CN.hpp.

Definition at line 133 of file CN.hpp.

Initial value:

this const vector defines the starting and ending EntityType for each dimension, e.g. TypeDimensionMap[2] returns a pair of EntityTypes bounding dimension 2.

Definition at line 138 of file CN.hpp.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines