Loading [MathJax]/extensions/tex2jax.js
Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
iMOAB.h File Reference
#include "moab/MOABConfig.h"
#include "moab/Types.hpp"
#include <stdio.h>
#include <assert.h>
#include <stdlib.h>
+ Include dependency graph for iMOAB.h:
+ This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Macros

#define iMOAB_AppID   int*
 
#define iMOAB_String   char*
 
#define iMOAB_GlobalID   int
 
#define iMOAB_LocalID   int
 
#define ErrCode   int
 
#define __PRETTY_FUNCTION__   __func__
 
#define CHK_MPI_ERR(ierr)
 
#define IMOAB_CHECKPOINTER(prmObj, position)
 
#define IMOAB_THROW_ERROR(message, retval)
 
#define IMOAB_ASSERT_RET(condition, message, retval)
 
#define IMOAB_ASSERT(condition, message)   IMOAB_ASSERT_RET( condition, message, moab::MB_UNHANDLED_OPTION )
 

Enumerations

enum  MOAB_TAG_TYPE {
  DENSE_INTEGER = 0 , DENSE_DOUBLE = 1 , DENSE_ENTITYHANDLE = 2 , SPARSE_INTEGER = 3 ,
  SPARSE_DOUBLE = 4 , SPARSE_ENTITYHANDLE = 5
}
 MOAB tag types can be: dense/sparse, and of int/double/EntityHandle types. They can also be defined on both elements and vertices of the mesh. More...
 
enum  MOAB_TAG_OWNER_TYPE { TAG_VERTEX = 0 , TAG_EDGE = 1 , TAG_FACE = 2 , TAG_ELEMENT = 3 }
 Define MOAB tag ownership information i.e., the entity where the tag is primarily defined This is typically used to query the tag data and to set it back. More...
 

Functions

ErrCode iMOAB_Initialize (int argc, iMOAB_String *argv)
 Initialize the iMOAB interface implementation. More...
 
ErrCode iMOAB_InitializeFortran (void)
 Initialize the iMOAB interface implementation from Fortran driver. More...
 
ErrCode iMOAB_Finalize (void)
 Finalize the iMOAB interface implementation. More...
 
ErrCode iMOAB_RegisterApplication (const iMOAB_String app_name, int *compid, iMOAB_AppID pid)
 Register application - Create a unique application ID and bootstrap interfaces for further queries. More...
 
ErrCode iMOAB_DeregisterApplication (iMOAB_AppID pid)
 De-Register application: delete mesh (set) associated with the application ID. More...
 
ErrCode iMOAB_RegisterApplicationFortran (const iMOAB_String app_name, int *compid, iMOAB_AppID pid)
 Register a Fortran-based application - Create a unique application ID and bootstrap interfaces for further queries. More...
 
ErrCode iMOAB_DeregisterApplicationFortran (iMOAB_AppID pid)
 De-Register the Fortran based application: delete mesh (set) associated with the application ID. More...
 
ErrCode iMOAB_ReadHeaderInfo (const iMOAB_String filename, int *num_global_vertices, int *num_global_elements, int *num_dimension, int *num_parts)
 It should be called on master task only, and information obtained could be broadcasted by the user. It is a fast lookup in the header of the file. More...
 
ErrCode iMOAB_LoadMesh (iMOAB_AppID pid, const iMOAB_String filename, const iMOAB_String read_options, int *num_ghost_layers)
 Load a MOAB mesh file in parallel and exchange ghost layers as requested. More...
 
ErrCode iMOAB_CreateVertices (iMOAB_AppID pid, int *coords_len, int *dim, double *coordinates)
 Create vertices for an app; it assumes no other vertices. More...
 
ErrCode iMOAB_CreateElements (iMOAB_AppID pid, int *num_elem, int *type, int *num_nodes_per_element, int *connectivity, int *block_ID)
 Create elements for an app; it assumes connectivity from local vertices, in order. More...
 
ErrCode iMOAB_ResolveSharedEntities (iMOAB_AppID pid, int *num_verts, int *marker)
 Resolve shared entities using global markers on shared vertices. More...
 
ErrCode iMOAB_DetermineGhostEntities (iMOAB_AppID pid, int *ghost_dim, int *num_ghost_layers, int *bridge_dim)
 Create the requested number of ghost layers for the parallel mesh. More...
 
ErrCode iMOAB_WriteMesh (iMOAB_AppID pid, const iMOAB_String filename, const iMOAB_String write_options)
 Write a MOAB mesh along with the solution tags to a file. More...
 
ErrCode iMOAB_WriteLocalMesh (iMOAB_AppID pid, iMOAB_String prefix)
 Write a local MOAB mesh copy. More...
 
ErrCode iMOAB_UpdateMeshInfo (iMOAB_AppID pid)
 Update local mesh data structure, from file information. More...
 
ErrCode iMOAB_GetMeshInfo (iMOAB_AppID pid, int *num_visible_vertices, int *num_visible_elements, int *num_visible_blocks, int *num_visible_surfaceBC, int *num_visible_vertexBC)
 Retrieve all the important mesh information related to vertices, elements, vertex and surface boundary conditions. More...
 
ErrCode iMOAB_GetVertexID (iMOAB_AppID pid, int *vertices_length, iMOAB_GlobalID *global_vertex_ID)
 Get the global vertex IDs for all locally visible (owned and shared/ghosted) vertices. More...
 
ErrCode iMOAB_GetVertexOwnership (iMOAB_AppID pid, int *vertices_length, int *visible_global_rank_ID)
 Get vertex ownership information. More...
 
ErrCode iMOAB_GetVisibleVerticesCoordinates (iMOAB_AppID pid, int *coords_length, double *coordinates)
 Get vertex coordinates for all local (owned and ghosted) vertices. More...
 
ErrCode iMOAB_GetBlockID (iMOAB_AppID pid, int *block_length, iMOAB_GlobalID *global_block_IDs)
 Get the global block IDs for all locally visible (owned and shared/ghosted) blocks. More...
 
ErrCode iMOAB_GetBlockInfo (iMOAB_AppID pid, iMOAB_GlobalID *global_block_ID, int *vertices_per_element, int *num_elements_in_block)
 Get the global block information and number of visible elements of belonging to a block (MATERIAL SET). More...
 
ErrCode iMOAB_GetVisibleElementsInfo (iMOAB_AppID pid, int *num_visible_elements, iMOAB_GlobalID *element_global_IDs, int *ranks, iMOAB_GlobalID *block_IDs)
 Get the visible elements information. More...
 
ErrCode iMOAB_GetBlockElementConnectivities (iMOAB_AppID pid, iMOAB_GlobalID *global_block_ID, int *connectivity_length, int *element_connectivity)
 Get the connectivities for elements contained within a certain block. More...
 
ErrCode iMOAB_GetElementConnectivity (iMOAB_AppID pid, iMOAB_LocalID *elem_index, int *connectivity_length, int *element_connectivity)
 Get the connectivity for one element only. More...
 
ErrCode iMOAB_GetElementOwnership (iMOAB_AppID pid, iMOAB_GlobalID *global_block_ID, int *num_elements_in_block, int *element_ownership)
 Get the element ownership within a certain block i.e., processor ID of the element owner. More...
 
ErrCode iMOAB_GetElementID (iMOAB_AppID pid, iMOAB_GlobalID *global_block_ID, int *num_elements_in_block, iMOAB_GlobalID *global_element_ID, iMOAB_LocalID *local_element_ID)
 Get the global IDs for all locally visible elements belonging to a particular block. More...
 
ErrCode iMOAB_GetPointerToSurfaceBC (iMOAB_AppID pid, int *surface_BC_length, iMOAB_LocalID *local_element_ID, int *reference_surface_ID, int *boundary_condition_value)
 Get the surface boundary condition information. More...
 
ErrCode iMOAB_GetPointerToVertexBC (iMOAB_AppID pid, int *vertex_BC_length, iMOAB_LocalID *local_vertex_ID, int *boundary_condition_value)
 Get the vertex boundary condition information. More...
 
ErrCode iMOAB_DefineTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_name, int *tag_type, int *components_per_entity, int *tag_index)
 Define a MOAB Tag corresponding to the application depending on requested types. More...
 
ErrCode iMOAB_SetIntTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_name, int *num_tag_storage_length, int *entity_type, int *tag_storage_data)
 Store the specified values in a MOAB integer Tag. More...
 
ErrCode iMOAB_GetIntTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_name, int *num_tag_storage_length, int *entity_type, int *tag_storage_data)
 Get the specified values in a MOAB integer Tag. More...
 
ErrCode iMOAB_SetDoubleTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_name, int *num_tag_storage_length, int *entity_type, double *tag_storage_data)
 Store the specified values in a MOAB double Tag. More...
 
ErrCode iMOAB_SetDoubleTagStorageWithGid (iMOAB_AppID pid, const iMOAB_String tag_storage_names, int *num_tag_storage_length, int *entity_type, double *tag_storage_data, int *globalIds)
 Store the specified values in a MOAB double Tag, for given ids. More...
 
ErrCode iMOAB_GetDoubleTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_name, int *num_tag_storage_length, int *entity_type, double *tag_storage_data)
 Retrieve the specified values in a MOAB double Tag. More...
 
ErrCode iMOAB_SynchronizeTags (iMOAB_AppID pid, int *num_tag, int *tag_indices, int *ent_type)
 Exchange tag values for the given tags across process boundaries. More...
 
ErrCode iMOAB_ReduceTagsMax (iMOAB_AppID pid, int *tag_index, int *ent_type)
 Perform global reductions with the processes in the current applications communicator. Specifically this routine performs reductions on the maximum value of the tag indicated by tag_index. More...
 
ErrCode iMOAB_GetNeighborElements (iMOAB_AppID pid, iMOAB_LocalID *local_index, int *num_adjacent_elements, iMOAB_LocalID *adjacent_element_IDs)
 Retrieve the adjacencies for the element entities. More...
 
ErrCode iMOAB_GetNeighborVertices (iMOAB_AppID pid, iMOAB_LocalID *local_vertex_ID, int *num_adjacent_vertices, iMOAB_LocalID *adjacent_vertex_IDs)
 Get the adjacencies for the vertex entities. More...
 
ErrCode iMOAB_SetGlobalInfo (iMOAB_AppID pid, int *num_global_verts, int *num_global_elems)
 Set global information for number of vertices and number of elements; It is usually available from the h5m file, or it can be computed with a MPI_Reduce. More...
 
ErrCode iMOAB_GetGlobalInfo (iMOAB_AppID pid, int *num_global_verts, int *num_global_elems)
 Get global information about number of vertices and number of elements. More...
 

Detailed Description

iMOAB: a language-agnostic, lightweight interface to MOAB.

Supports usage from C/C++, Fortran (77/90/2003), Python.

Remarks
1) All data in the interface are exposed via POD-types.
2) Pass everything by reference, so we do not have to use VAL()
3) Arrays are allocated by the user code. No concerns about de-allocation of the data will be taken up by the interface.
4) Always pass the pointer to the start of array along with the total allocated size for the array.
5) Return the filled array requested by user along with optionally the actual length of the array that was filled. (for typical cases, should be the allocated length)

Definition in file iMOAB.h.

Macro Definition Documentation

◆ __PRETTY_FUNCTION__

#define __PRETTY_FUNCTION__   __func__

Definition at line 65 of file iMOAB.h.

◆ CHK_MPI_ERR

#define CHK_MPI_ERR (   ierr)
Value:
{ \ if( 0 != ( ierr ) ) return moab::MB_FAILURE; \ }

Definition at line 95 of file iMOAB.h.

◆ ErrCode

#define ErrCode   int

Definition at line 58 of file iMOAB.h.

◆ iMOAB_AppID

#define iMOAB_AppID   int*

\Notes 1) Fortran MPI_Comm won't work. Take an integer argument and use MPI_F2C calls to get the C-Comm object 2) ReadHeaderInfo - Does it need the pid ? 3) Reuse the comm object from the registration for both load and write operations. Do not take comm objects again to avoid confusion and retain consistency. 4) Decipher the global comm object and the subset partioning for each application based on the comm object 5) GetMeshInfo - return separately the owned and ghosted vertices/elements – not together in visible_** but rather owned_** and ghosted_**. Make these arrays of size 2. 6) Should we sort the vertices after ghosting, such that the locally owned is first, and the ghosted appended next. 7) RCM only for the owned part of the mesh – do not screw with the ghosted layers 8) GetBlockID - identical to GetVertexID – return the global numbering for block 9) GetVertexID – remember that the order of vertices returned have an implicit numbering embedded in it. DO NOT CHANGE THIS ORDERING... 10) GetBlockInfo takes global Block ID; Remove blockname unless there is a separate use case for it.. 11) GetElementConnectivity - clarify whether we return global or local vertex numbering. Preferably local numbering else lot of deciphering for global.

Definition at line 51 of file iMOAB.h.

◆ IMOAB_ASSERT

#define IMOAB_ASSERT (   condition,
  message 
)    IMOAB_ASSERT_RET( condition, message, moab::MB_UNHANDLED_OPTION )

Definition at line 126 of file iMOAB.h.

◆ IMOAB_ASSERT_RET

#define IMOAB_ASSERT_RET (   condition,
  message,
  retval 
)
Value:
{ \ if( !( condition ) ) \ { \ printf( "iMOAB Error: %s.\n Failed condition: %s\n Location: %s in %s:%d.\n", message, #condition, \ __PRETTY_FUNCTION__, __FILE__, __LINE__ ); \ return retval; \ } \ }

Definition at line 116 of file iMOAB.h.

◆ IMOAB_CHECKPOINTER

#define IMOAB_CHECKPOINTER (   prmObj,
  position 
)
Value:
{ \ if( prmObj == nullptr ) \ { \ printf( "InputParamError at %d: '%s' is invalid and null.\n", position, #prmObj ); \ return moab::MB_UNHANDLED_OPTION; \ } \ }

Definition at line 100 of file iMOAB.h.

◆ iMOAB_GlobalID

#define iMOAB_GlobalID   int

Definition at line 53 of file iMOAB.h.

◆ iMOAB_LocalID

#define iMOAB_LocalID   int

Definition at line 54 of file iMOAB.h.

◆ iMOAB_String

#define iMOAB_String   char*

Definition at line 52 of file iMOAB.h.

◆ IMOAB_THROW_ERROR

#define IMOAB_THROW_ERROR (   message,
  retval 
)
Value:
{ \ printf( "iMOAB Error: %s.\n Location: %s in %s:%d.\n", message, __PRETTY_FUNCTION__, __FILE__, __LINE__ ); \ return retval; \ }

Definition at line 109 of file iMOAB.h.

Enumeration Type Documentation

◆ MOAB_TAG_OWNER_TYPE

Define MOAB tag ownership information i.e., the entity where the tag is primarily defined This is typically used to query the tag data and to set it back.

Enumerator
TAG_VERTEX 
TAG_EDGE 
TAG_FACE 
TAG_ELEMENT 

Definition at line 87 of file iMOAB.h.

88 { 89  TAG_VERTEX = 0, 90  TAG_EDGE = 1, 91  TAG_FACE = 2, 92  TAG_ELEMENT = 3 93 };

◆ MOAB_TAG_TYPE

MOAB tag types can be: dense/sparse, and of int/double/EntityHandle types. They can also be defined on both elements and vertices of the mesh.

Enumerator
DENSE_INTEGER 
DENSE_DOUBLE 
DENSE_ENTITYHANDLE 
SPARSE_INTEGER 
SPARSE_DOUBLE 
SPARSE_ENTITYHANDLE 

Definition at line 73 of file iMOAB.h.

74 { 75  DENSE_INTEGER = 0, 76  DENSE_DOUBLE = 1, 77  DENSE_ENTITYHANDLE = 2, 78  SPARSE_INTEGER = 3, 79  SPARSE_DOUBLE = 4, 80  SPARSE_ENTITYHANDLE = 5 81 };

Function Documentation

◆ iMOAB_CreateElements()

ErrCode iMOAB_CreateElements ( iMOAB_AppID  pid,
int *  num_elem,
int *  type,
int *  num_nodes_per_element,
int *  connectivity,
int *  block_ID 
)

Create elements for an app; it assumes connectivity from local vertices, in order.

Note
Operations: Not collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]num_elem(int*) Number of elements.
[in]type(int*) Type of element (moab type).
[in]num_nodes_per_element(int*) Number of nodes per element.
[in]connectivity(int *) Connectivity array, with respect to vertices (assumed contiguous).
[in]block_ID(int *) The local block identifier, which will now contain the elements.
Returns
ErrCode The error code indicating success or failure.

add the new ents to the clock set

Definition at line 2243 of file iMOAB.cpp.

2249 { 2250  // Create elements 2251  appData& data = context.appDatas[*pid]; 2252  2253  ReadUtilIface* read_iface; 2254  ErrorCode rval = context.MBI->query_interface( read_iface );MB_CHK_ERR( rval ); 2255  2256  EntityType mbtype = (EntityType)( *type ); 2257  EntityHandle actual_start_handle; 2258  EntityHandle* array = NULL; 2259  rval = read_iface->get_element_connect( *num_elem, *num_nodes_per_element, mbtype, 1, actual_start_handle, array );MB_CHK_ERR( rval ); 2260  2261  // fill up with actual connectivity from input; assume the vertices are in order, and start 2262  // vertex is the first in the current data vertex range 2263  EntityHandle firstVertex = data.local_verts[0]; 2264  2265  for( int j = 0; j < *num_elem * ( *num_nodes_per_element ); j++ ) 2266  { 2267  array[j] = connectivity[j] + firstVertex - 1; 2268  } // assumes connectivity uses 1 based array (from fortran, mostly) 2269  2270  Range new_elems( actual_start_handle, actual_start_handle + *num_elem - 1 ); 2271  2272  rval = context.MBI->add_entities( data.file_set, new_elems );MB_CHK_ERR( rval ); 2273  2274  data.primary_elems.merge( new_elems ); 2275  2276  // add to adjacency 2277  rval = read_iface->update_adjacencies( actual_start_handle, *num_elem, *num_nodes_per_element, array );MB_CHK_ERR( rval ); 2278  // organize all new elements in block, with the given block ID; if the block set is not 2279  // existing, create a new mesh set; 2280  Range sets; 2281  int set_no = *block_ID; 2282  const void* setno_ptr = &set_no; 2283  rval = context.MBI->get_entities_by_type_and_tag( data.file_set, MBENTITYSET, &context.material_tag, &setno_ptr, 1, 2284  sets ); 2285  EntityHandle block_set; 2286  2287  if( MB_FAILURE == rval || sets.empty() ) 2288  { 2289  // create a new set, with this block ID 2290  rval = context.MBI->create_meshset( MESHSET_SET, block_set );MB_CHK_ERR( rval ); 2291  2292  rval = context.MBI->tag_set_data( context.material_tag, &block_set, 1, &set_no );MB_CHK_ERR( rval ); 2293  2294  // add the material set to file set 2295  rval = context.MBI->add_entities( data.file_set, &block_set, 1 );MB_CHK_ERR( rval ); 2296  } 2297  else 2298  { 2299  block_set = sets[0]; 2300  } // first set is the one we want 2301  2302  /// add the new ents to the clock set 2303  rval = context.MBI->add_entities( block_set, new_elems );MB_CHK_ERR( rval ); 2304  2305  return moab::MB_SUCCESS; 2306 }

References moab::Interface::add_entities(), GlobalContext::appDatas, context, moab::Interface::create_meshset(), moab::Range::empty(), ErrorCode, appData::file_set, moab::ReadUtilIface::get_element_connect(), moab::Interface::get_entities_by_type_and_tag(), appData::local_verts, GlobalContext::material_tag, MB_CHK_ERR, MB_SUCCESS, MBENTITYSET, GlobalContext::MBI, moab::Range::merge(), MESHSET_SET, appData::primary_elems, moab::Interface::query_interface(), moab::Interface::tag_set_data(), and moab::ReadUtilIface::update_adjacencies().

◆ iMOAB_CreateVertices()

ErrCode iMOAB_CreateVertices ( iMOAB_AppID  pid,
int *  coords_len,
int *  dim,
double *  coordinates 
)

Create vertices for an app; it assumes no other vertices.

Note
Operations: Not collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]coords_len(int*) Size of the coordinates array (nverts * dim).
[in]dim(int*) Dimension of the vertex coordinates (usually 3).
[in]coordinates(double*) Coordinates of all vertices, interleaved.
Returns
ErrCode The error code indicating success or failure.

Definition at line 2222 of file iMOAB.cpp.

2223 { 2224  ErrorCode rval; 2225  appData& data = context.appDatas[*pid]; 2226  2227  if( !data.local_verts.empty() ) // we should have no vertices in the app 2228  { 2229  return moab::MB_FAILURE; 2230  } 2231  2232  int nverts = *coords_len / *dim; 2233  2234  rval = context.MBI->create_vertices( coordinates, nverts, data.local_verts );MB_CHK_ERR( rval ); 2235  2236  rval = context.MBI->add_entities( data.file_set, data.local_verts );MB_CHK_ERR( rval ); 2237  2238  // also add the vertices to the all_verts range 2239  data.all_verts.merge( data.local_verts ); 2240  return moab::MB_SUCCESS; 2241 }

References moab::Interface::add_entities(), appData::all_verts, GlobalContext::appDatas, context, moab::Interface::create_vertices(), dim, moab::Range::empty(), ErrorCode, appData::file_set, appData::local_verts, MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, and moab::Range::merge().

◆ iMOAB_DefineTagStorage()

ErrCode iMOAB_DefineTagStorage ( iMOAB_AppID  pid,
const iMOAB_String  tag_storage_name,
int *  tag_type,
int *  components_per_entity,
int *  tag_index 
)

Define a MOAB Tag corresponding to the application depending on requested types.

Note
In MOAB, for most solution vectors, we only need to create a "Dense", "Double" Tag. A sparse tag can be created too. If the tag is already existing in the file, it will not be created. If it is a new tag, memory will be allocated when setting the values. Default values are 0 for for integer tags, 0.0 for double tags, 0 for entity handle tags.

Operations: Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]tag_storage_name(iMOAB_String) The tag name to store/retrieve the data in MOAB.
[in]tag_type(int*) The type of MOAB tag (Dense/Sparse, Double/Int/EntityHandle), enum MOAB_TAG_TYPE.
[in]components_per_entity(int*) The total size of vector dimension per entity for the tag (e.g., number of doubles per entity).
[out]tag_index(int*) The tag index which can be used as identifier in synchronize methods.
Returns
ErrCode The error code indicating success or failure.

Definition at line 1471 of file iMOAB.cpp.

1476 { 1477  // see if the tag is already existing, and if yes, check the type, length 1478  if( *tag_type < 0 || *tag_type > 5 ) 1479  { 1480  return moab::MB_FAILURE; 1481  } // we have 6 types of tags supported so far 1482  1483  DataType tagDataType; 1484  TagType tagType; 1485  void* defaultVal = NULL; 1486  int* defInt = new int[*components_per_entity]; 1487  double* defDouble = new double[*components_per_entity]; 1488  EntityHandle* defHandle = new EntityHandle[*components_per_entity]; 1489  1490  for( int i = 0; i < *components_per_entity; i++ ) 1491  { 1492  defInt[i] = 0; 1493  defDouble[i] = -1e+10; 1494  defHandle[i] = (EntityHandle)0; 1495  } 1496  1497  switch( *tag_type ) 1498  { 1499  case 0: 1500  tagDataType = MB_TYPE_INTEGER; 1501  tagType = MB_TAG_DENSE; 1502  defaultVal = defInt; 1503  break; 1504  1505  case 1: 1506  tagDataType = MB_TYPE_DOUBLE; 1507  tagType = MB_TAG_DENSE; 1508  defaultVal = defDouble; 1509  break; 1510  1511  case 2: 1512  tagDataType = MB_TYPE_HANDLE; 1513  tagType = MB_TAG_DENSE; 1514  defaultVal = defHandle; 1515  break; 1516  1517  case 3: 1518  tagDataType = MB_TYPE_INTEGER; 1519  tagType = MB_TAG_SPARSE; 1520  defaultVal = defInt; 1521  break; 1522  1523  case 4: 1524  tagDataType = MB_TYPE_DOUBLE; 1525  tagType = MB_TAG_SPARSE; 1526  defaultVal = defDouble; 1527  break; 1528  1529  case 5: 1530  tagDataType = MB_TYPE_HANDLE; 1531  tagType = MB_TAG_SPARSE; 1532  defaultVal = defHandle; 1533  break; 1534  1535  default: { 1536  delete[] defInt; 1537  delete[] defDouble; 1538  delete[] defHandle; 1539  return moab::MB_FAILURE; 1540  } // error 1541  } 1542  1543  Tag tagHandle; 1544  // split storage names if separated list 1545  1546  std::string tag_name( tag_storage_name ); 1547  1548  // first separate the names of the tags 1549  // we assume that there are separators ":" between the tag names 1550  std::vector< std::string > tagNames; 1551  std::string separator( ":" ); 1552  split_tag_names( tag_name, separator, tagNames ); 1553  1554  ErrorCode rval = moab::MB_SUCCESS; // assume success already :) 1555  appData& data = context.appDatas[*pid]; 1556  int already_defined_tags = 0; 1557  1558  for( size_t i = 0; i < tagNames.size(); i++ ) 1559  { 1560  rval = context.MBI->tag_get_handle( tagNames[i].c_str(), *components_per_entity, tagDataType, tagHandle, 1561  tagType | MB_TAG_EXCL | MB_TAG_CREAT, defaultVal ); 1562  1563  if( MB_ALREADY_ALLOCATED == rval ) 1564  { 1565  std::map< std::string, Tag >& mTags = data.tagMap; 1566  std::map< std::string, Tag >::iterator mit = mTags.find( tagNames[i].c_str() ); 1567  1568  if( mit == mTags.end() ) 1569  { 1570  // add it to the map 1571  mTags[tagNames[i]] = tagHandle; 1572  // push it to the list of tags, too 1573  *tag_index = (int)data.tagList.size(); 1574  data.tagList.push_back( tagHandle ); 1575  } 1576  rval = MB_SUCCESS; 1577  already_defined_tags++; 1578  } 1579  else if( MB_SUCCESS == rval ) 1580  { 1581  data.tagMap[tagNames[i]] = tagHandle; 1582  *tag_index = (int)data.tagList.size(); 1583  data.tagList.push_back( tagHandle ); 1584  } 1585  else 1586  { 1587  rval = moab::MB_FAILURE; // some tags were not created 1588  } 1589  } 1590  // we don't need default values anymore, avoid leaks 1591  int rankHere = 0; 1592 #ifdef MOAB_HAVE_MPI 1593  ParallelComm* pco = context.pcomms[*pid]; 1594  rankHere = pco->rank(); 1595 #endif 1596  if( !rankHere && already_defined_tags ) 1597  std::cout << " application with ID: " << *pid << " global id: " << data.global_id << " name: " << data.name 1598  << " has " << already_defined_tags << " already defined tags out of " << tagNames.size() 1599  << " tags \n"; 1600  delete[] defInt; 1601  delete[] defDouble; 1602  delete[] defHandle; 1603  return rval; 1604 }

References GlobalContext::appDatas, context, ErrorCode, appData::global_id, MB_ALREADY_ALLOCATED, MB_SUCCESS, MB_TAG_CREAT, MB_TAG_DENSE, MB_TAG_EXCL, MB_TAG_SPARSE, MB_TYPE_DOUBLE, MB_TYPE_HANDLE, MB_TYPE_INTEGER, GlobalContext::MBI, appData::name, GlobalContext::pcomms, moab::ParallelComm::rank(), split_tag_names(), moab::Interface::tag_get_handle(), appData::tagList, appData::tagMap, and TagType.

Referenced by iMOAB_ReceiveMesh().

◆ iMOAB_DeregisterApplication()

ErrCode iMOAB_DeregisterApplication ( iMOAB_AppID  pid)

De-Register application: delete mesh (set) associated with the application ID.

Note
The associated communicator will be released, and all associated mesh entities and sets will be deleted from the mesh data structure. Associated tag storage data will be freed too.

Operations: Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID
Returns
ErrCode The error code indicating success or failure.

Definition at line 335 of file iMOAB.cpp.

336 { 337  // the file set , parallel comm are all in vectors indexed by *pid 338  // assume we did not delete anything yet 339  // *pid will not be reused if we register another application 340  appData& data = context.appDatas[*pid]; 341  int rankHere = 0; 342 #ifdef MOAB_HAVE_MPI 343  ParallelComm* pco = context.pcomms[*pid]; 344  rankHere = pco->rank(); 345 #endif 346  if( !rankHere ) 347  std::cout << " application with ID: " << *pid << " global id: " << data.global_id << " name: " << data.name 348  << " is de-registered now \n"; 349  350  EntityHandle fileSet = data.file_set; 351  // get all entities part of the file set 352  Range fileents; 353  ErrorCode rval = context.MBI->get_entities_by_handle( fileSet, fileents, /*recursive */ true );MB_CHK_ERR( rval ); 354  355  fileents.insert( fileSet ); 356  357  rval = context.MBI->get_entities_by_type( fileSet, MBENTITYSET, fileents );MB_CHK_ERR( rval ); // append all mesh sets 358  359 #ifdef MOAB_HAVE_TEMPESTREMAP 360  if( data.tempestData.remapper ) delete data.tempestData.remapper; 361  if( data.tempestData.weightMaps.size() ) data.tempestData.weightMaps.clear(); 362 #endif 363  364 #ifdef MOAB_HAVE_MPI 365  366  // we could get the pco also with 367  // ParallelComm * pcomm = ParallelComm::get_pcomm(context.MBI, *pid); 368  369  std::map< int, ParCommGraph* >& pargs = data.pgraph; 370  371  // free the parallel comm graphs associated with this app 372  for( std::map< int, ParCommGraph* >::iterator mt = pargs.begin(); mt != pargs.end(); mt++ ) 373  { 374  ParCommGraph* pgr = mt->second; 375  delete pgr; 376  pgr = NULL; 377  } 378  if( pco ) delete pco; 379 #endif 380  381  // delete first all except vertices 382  Range vertices = fileents.subset_by_type( MBVERTEX ); 383  Range noverts = subtract( fileents, vertices ); 384  385  rval = context.MBI->delete_entities( noverts );MB_CHK_ERR( rval ); 386  // now retrieve connected elements that still exist (maybe in other sets, pids?) 387  Range adj_ents_left; 388  rval = context.MBI->get_adjacencies( vertices, 1, false, adj_ents_left, Interface::UNION );MB_CHK_ERR( rval ); 389  rval = context.MBI->get_adjacencies( vertices, 2, false, adj_ents_left, Interface::UNION );MB_CHK_ERR( rval ); 390  rval = context.MBI->get_adjacencies( vertices, 3, false, adj_ents_left, Interface::UNION );MB_CHK_ERR( rval ); 391  392  if( !adj_ents_left.empty() ) 393  { 394  Range conn_verts; 395  rval = context.MBI->get_connectivity( adj_ents_left, conn_verts );MB_CHK_ERR( rval ); 396  vertices = subtract( vertices, conn_verts ); 397  } 398  399  rval = context.MBI->delete_entities( vertices );MB_CHK_ERR( rval ); 400  401  std::map< std::string, int >::iterator mit; 402  403  for( mit = context.appIdMap.begin(); mit != context.appIdMap.end(); mit++ ) 404  { 405  int pidx = mit->second; 406  407  if( *pid == pidx ) 408  { 409  break; 410  } 411  } 412  413  context.appIdMap.erase( mit ); 414  std::map< int, int >::iterator mit1; 415  416  for( mit1 = context.appIdCompMap.begin(); mit1 != context.appIdCompMap.end(); mit1++ ) 417  { 418  int pidx = mit1->second; 419  420  if( *pid == pidx ) 421  { 422  break; 423  } 424  } 425  426  context.appIdCompMap.erase( mit1 ); 427  428  context.unused_pid--; // we have to go backwards always TODO 429  context.appDatas.pop_back(); 430 #ifdef MOAB_HAVE_MPI 431  context.pcomms.pop_back(); 432 #endif 433  return moab::MB_SUCCESS; 434 }

References GlobalContext::appDatas, GlobalContext::appIdCompMap, GlobalContext::appIdMap, context, moab::Interface::delete_entities(), moab::Range::empty(), ErrorCode, appData::file_set, moab::Interface::get_adjacencies(), moab::Interface::get_connectivity(), moab::Interface::get_entities_by_handle(), moab::Interface::get_entities_by_type(), appData::global_id, moab::Range::insert(), MB_CHK_ERR, MB_SUCCESS, MBENTITYSET, GlobalContext::MBI, MBVERTEX, appData::name, GlobalContext::pcomms, appData::pgraph, moab::ParallelComm::rank(), moab::Range::subset_by_type(), moab::subtract(), moab::Interface::UNION, and GlobalContext::unused_pid.

Referenced by iMOAB_DeregisterApplicationFortran().

◆ iMOAB_DeregisterApplicationFortran()

ErrCode iMOAB_DeregisterApplicationFortran ( iMOAB_AppID  pid)

De-Register the Fortran based application: delete mesh (set) associated with the application ID.

Note
The associated communicator will be released, and all associated mesh entities and sets will be deleted from the mesh data structure. Associated tag storage data will be freed too.

Operations: Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID
Returns
ErrCode The error code indicating success or failure.

Definition at line 436 of file iMOAB.cpp.

437 { 438  // release any Fortran specific allocations here before we pass it on 439  context.appDatas[*pid].is_fortran = false; 440  441  // release all datastructure allocations 442  return iMOAB_DeregisterApplication( pid ); 443 }

References GlobalContext::appDatas, context, and iMOAB_DeregisterApplication().

◆ iMOAB_DetermineGhostEntities()

ErrCode iMOAB_DetermineGhostEntities ( iMOAB_AppID  pid,
int *  ghost_dim,
int *  num_ghost_layers,
int *  bridge_dim 
)

Create the requested number of ghost layers for the parallel mesh.

Note
The routine assumes that the shared entities were resolved successfully, and that the mesh is properly distributed on separate tasks.

Operations: Collective.

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]ghost_dim(int*) Desired ghost dimension (2 or 3, most of the time).
[in]num_ghost_layers(int*) Number of ghost layers requested.
[in]bridge_dim(int*) Bridge dimension (0 for vertices, 1 for edges, 2 for faces).
Returns
ErrCode The error code indicating success or failure.

Definition at line 2385 of file iMOAB.cpp.

2386 { 2387  ErrorCode rval; 2388  2389  // verify we have valid ghost layers input specified. If invalid, exit quick. 2390  if( num_ghost_layers && *num_ghost_layers <= 0 ) 2391  { 2392  return moab::MB_SUCCESS; 2393  } // nothing to do 2394  2395  appData& data = context.appDatas[*pid]; 2396  ParallelComm* pco = context.pcomms[*pid]; 2397  2398  // maybe we should be passing this too; most of the time we do not need additional ents collective call 2399  constexpr int addl_ents = 0; 2400  rval = pco->exchange_ghost_cells( *ghost_dim, *bridge_dim, 1, // get only one layer of ghost entities 2401  addl_ents, true, true, &data.file_set );MB_CHK_ERR( rval ); 2402  for( int i = 2; i <= *num_ghost_layers; i++ ) 2403  { 2404  rval = pco->correct_thin_ghost_layers();MB_CHK_ERR( rval ); // correct for thin layers 2405  rval = pco->exchange_ghost_cells( *ghost_dim, *bridge_dim, i, // iteratively get one extra layer 2406  addl_ents, true, true, &data.file_set );MB_CHK_ERR( rval ); 2407  } 2408  2409  // Update ghost layer information 2410  data.num_ghost_layers = *num_ghost_layers; 2411  2412  // now re-establish all mesh info; will reconstruct mesh info, based solely on what is in the file set 2413  return iMOAB_UpdateMeshInfo( pid ); 2414 }

References GlobalContext::appDatas, context, moab::ParallelComm::correct_thin_ghost_layers(), ErrorCode, moab::ParallelComm::exchange_ghost_cells(), appData::file_set, iMOAB_UpdateMeshInfo(), MB_CHK_ERR, MB_SUCCESS, appData::num_ghost_layers, and GlobalContext::pcomms.

◆ iMOAB_Finalize()

ErrCode iMOAB_Finalize ( void  )

Finalize the iMOAB interface implementation.

It will delete the internally reference counted MOAB instance, if the reference count reaches 0.

Operations: Collective

Returns
ErrCode The error code indicating success or failure.

Definition at line 189 of file iMOAB.cpp.

190 { 191  context.refCountMB--; 192  193  if( 0 == context.refCountMB ) 194  { 195  delete context.MBI; 196  } 197  198  return MB_SUCCESS; 199 }

References context, MB_SUCCESS, GlobalContext::MBI, and GlobalContext::refCountMB.

◆ iMOAB_GetBlockElementConnectivities()

ErrCode iMOAB_GetBlockElementConnectivities ( iMOAB_AppID  pid,
iMOAB_GlobalID global_block_ID,
int *  connectivity_length,
int *  element_connectivity 
)

Get the connectivities for elements contained within a certain block.

Note
input is the block ID. Should we change to visible block local index?

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]global_block_ID(iMOAB_GlobalID*) The global block ID of the set being queried.
[in]connectivity_length(int*) The allocated size of array. (typical size := vertices_per_element*num_elements_in_block)
[out]element_connectivity(int*) The connectivity array to store element ordering in MOAB canonical numbering scheme (array allocated by user); array contains vertex indices in the local numbering order for vertices elements are in the same order as provided by GetElementOwnership and GetElementID.
Returns
ErrCode The error code indicating success or failure.

Definition at line 1163 of file iMOAB.cpp.

1167 { 1168  assert( global_block_ID ); // ensure global block ID argument is not null 1169  assert( connectivity_length ); // ensure connectivity length argument is not null 1170  1171  appData& data = context.appDatas[*pid]; 1172  std::map< int, int >& matMap = data.matIndex; 1173  std::map< int, int >::iterator it = matMap.find( *global_block_ID ); 1174  1175  if( it == matMap.end() ) 1176  { 1177  return moab::MB_FAILURE; 1178  } // error in finding block with id 1179  1180  int blockIndex = matMap[*global_block_ID]; 1181  EntityHandle matMeshSet = data.mat_sets[blockIndex]; 1182  std::vector< EntityHandle > elems; 1183  1184  ErrorCode rval = context.MBI->get_entities_by_handle( matMeshSet, elems );MB_CHK_ERR( rval ); 1185  1186  if( elems.empty() ) 1187  { 1188  return moab::MB_FAILURE; 1189  } 1190  1191  std::vector< EntityHandle > vconnect; 1192  rval = context.MBI->get_connectivity( &elems[0], elems.size(), vconnect );MB_CHK_ERR( rval ); 1193  1194  if( *connectivity_length != (int)vconnect.size() ) 1195  { 1196  return moab::MB_FAILURE; 1197  } // mismatched sizes 1198  1199  for( int i = 0; i < *connectivity_length; i++ ) 1200  { 1201  int inx = data.all_verts.index( vconnect[i] ); 1202  1203  if( -1 == inx ) 1204  { 1205  return moab::MB_FAILURE; 1206  } // error, vertex not in local range 1207  1208  element_connectivity[i] = inx; 1209  } 1210  1211  return moab::MB_SUCCESS; 1212 }

References appData::all_verts, GlobalContext::appDatas, context, ErrorCode, moab::Interface::get_connectivity(), moab::Interface::get_entities_by_handle(), moab::Range::index(), appData::mat_sets, appData::matIndex, MB_CHK_ERR, MB_SUCCESS, and GlobalContext::MBI.

◆ iMOAB_GetBlockID()

ErrCode iMOAB_GetBlockID ( iMOAB_AppID  pid,
int *  block_length,
iMOAB_GlobalID global_block_IDs 
)

Get the global block IDs for all locally visible (owned and shared/ghosted) blocks.

Note
Block IDs are corresponding to MATERIAL_SET tags for material sets. Usually the block ID is exported from Cubit as a unique integer value. First blocks are local, and next blocks are fully ghosted. First blocks have at least one owned cell/element, ghost blocks have only ghost cells. Internally, a block corresponds to a mesh set with a MATERIAL_SET tag value equal to the block ID.

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]block_length(int*) The allocated size of array (typical size := num_visible_blocks).
[out]global_block_IDs(iMOAB_GlobalID*) The global IDs for all locally visible blocks (array allocated by user).
Returns
ErrCode The error code indicating success or failure.

Definition at line 1040 of file iMOAB.cpp.

1041 { 1042  // local id blocks? they are counted from 0 to number of visible blocks ... 1043  // will actually return material set tag value for global 1044  Range& matSets = context.appDatas[*pid].mat_sets; 1045  1046  if( *block_length != (int)matSets.size() ) 1047  { 1048  return moab::MB_FAILURE; 1049  } 1050  1051  // return material set tag gtags[0 is material set tag 1052  ErrorCode rval = context.MBI->tag_get_data( context.material_tag, matSets, global_block_IDs );MB_CHK_ERR( rval ); 1053  1054  // populate map with index 1055  std::map< int, int >& matIdx = context.appDatas[*pid].matIndex; 1056  for( unsigned i = 0; i < matSets.size(); i++ ) 1057  { 1058  matIdx[global_block_IDs[i]] = i; 1059  } 1060  1061  return moab::MB_SUCCESS; 1062 }

References GlobalContext::appDatas, context, ErrorCode, GlobalContext::material_tag, MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, moab::Range::size(), and moab::Interface::tag_get_data().

◆ iMOAB_GetBlockInfo()

ErrCode iMOAB_GetBlockInfo ( iMOAB_AppID  pid,
iMOAB_GlobalID global_block_ID,
int *  vertices_per_element,
int *  num_elements_in_block 
)

Get the global block information and number of visible elements of belonging to a block (MATERIAL SET).

Note
A block has to be homogeneous, it can contain elements of a single type.

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]global_block_ID(iMOAB_GlobalID) The global block ID of the set to be queried.
[out]vertices_per_element(int*) The number of vertices per element.
[out]num_elements_in_block(int*) The number of elements in block.
Returns
ErrCode The error code indicating success or failure.

Definition at line 1064 of file iMOAB.cpp.

1068 { 1069  assert( global_block_ID ); 1070  1071  std::map< int, int >& matMap = context.appDatas[*pid].matIndex; 1072  std::map< int, int >::iterator it = matMap.find( *global_block_ID ); 1073  1074  if( it == matMap.end() ) 1075  { 1076  return moab::MB_FAILURE; 1077  } // error in finding block with id 1078  1079  int blockIndex = matMap[*global_block_ID]; 1080  EntityHandle matMeshSet = context.appDatas[*pid].mat_sets[blockIndex]; 1081  Range blo_elems; 1082  ErrorCode rval = context.MBI->get_entities_by_handle( matMeshSet, blo_elems ); 1083  1084  if( MB_SUCCESS != rval || blo_elems.empty() ) 1085  { 1086  return moab::MB_FAILURE; 1087  } 1088  1089  EntityType type = context.MBI->type_from_handle( blo_elems[0] ); 1090  1091  if( !blo_elems.all_of_type( type ) ) 1092  { 1093  return moab::MB_FAILURE; 1094  } // not all of same type 1095  1096  const EntityHandle* conn = NULL; 1097  int num_verts = 0; 1098  rval = context.MBI->get_connectivity( blo_elems[0], conn, num_verts );MB_CHK_ERR( rval ); 1099  1100  *vertices_per_element = num_verts; 1101  *num_elements_in_block = (int)blo_elems.size(); 1102  1103  return moab::MB_SUCCESS; 1104 }

References moab::Range::all_of_type(), GlobalContext::appDatas, context, moab::Range::empty(), ErrorCode, moab::Interface::get_connectivity(), moab::Interface::get_entities_by_handle(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, moab::Range::size(), and moab::Interface::type_from_handle().

◆ iMOAB_GetDoubleTagStorage()

ErrCode iMOAB_GetDoubleTagStorage ( iMOAB_AppID  pid,
const iMOAB_String  tag_storage_name,
int *  num_tag_storage_length,
int *  entity_type,
double *  tag_storage_data 
)

Retrieve the specified values in a MOAB double Tag.

Note
Operations: Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]tag_storage_name(iMOAB_String) The tag names to retrieve the data in MOAB.
[in]num_tag_storage_length(int) The size of total tag storage data (e.g., num_visible_vertices*components_per_entity*num_tags or num_visible_elements*components_per_entity*num_tags)
[in]entity_type(int*) Type=0 for vertices, and Type=1 for primary elements.
[out]tag_storage_data(double*) The array data of type double to replace the internal tag memory; The data is assumed to be contiguous over the local set of visible entities (either vertices or elements). unrolled by tags
Returns
ErrCode The error code indicating success or failure.

Definition at line 2039 of file iMOAB.cpp.

2044 { 2045  ErrorCode rval; 2046  // exactly the same code, except tag type check 2047  std::string tag_names( tag_storage_names ); 2048  // exactly the same code as for int tag :) maybe should check the type of tag too 2049  std::vector< std::string > tagNames; 2050  std::vector< Tag > tagHandles; 2051  std::string separator( ":" ); 2052  split_tag_names( tag_names, separator, tagNames ); 2053  2054  appData& data = context.appDatas[*pid]; 2055  2056  // set it on a subset of entities, based on type and length 2057  Range* ents_to_get = NULL; 2058  2059  if( *ent_type == 0 ) // vertices 2060  { 2061  ents_to_get = &data.all_verts; 2062  } 2063  else if( *ent_type == 1 ) 2064  { 2065  ents_to_get = &data.primary_elems; 2066  } 2067  int nents_to_get = (int)ents_to_get->size(); 2068  int position = 0; 2069  for( size_t i = 0; i < tagNames.size(); i++ ) 2070  { 2071  if( data.tagMap.find( tagNames[i] ) == data.tagMap.end() ) 2072  { 2073  return moab::MB_FAILURE; 2074  } // tag not defined 2075  2076  Tag tag = data.tagMap[tagNames[i]]; 2077  2078  int tagLength = 0; 2079  rval = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval ); 2080  2081  DataType dtype; 2082  rval = context.MBI->tag_get_data_type( tag, dtype );MB_CHK_ERR( rval ); 2083  2084  if( dtype != MB_TYPE_DOUBLE ) 2085  { 2086  return moab::MB_FAILURE; 2087  } 2088  2089  if( position + nents_to_get * tagLength > *num_tag_storage_length ) 2090  return moab::MB_FAILURE; // too many entity values to get 2091  2092  rval = context.MBI->tag_get_data( tag, *ents_to_get, &tag_storage_data[position] );MB_CHK_ERR( rval ); 2093  position = position + nents_to_get * tagLength; 2094  } 2095  2096  return moab::MB_SUCCESS; // no error 2097 }

References appData::all_verts, GlobalContext::appDatas, context, ErrorCode, MB_CHK_ERR, MB_SUCCESS, MB_TYPE_DOUBLE, GlobalContext::MBI, appData::primary_elems, moab::Range::size(), split_tag_names(), moab::Interface::tag_get_data(), moab::Interface::tag_get_data_type(), moab::Interface::tag_get_length(), and appData::tagMap.

◆ iMOAB_GetElementConnectivity()

ErrCode iMOAB_GetElementConnectivity ( iMOAB_AppID  pid,
iMOAB_LocalID elem_index,
int *  connectivity_length,
int *  element_connectivity 
)

Get the connectivity for one element only.

Note
This is a convenience method and in general should not be used unless necessary. You should use colaesced calls for multiple elements for efficiency.

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]elem_index(iMOAB_LocalID *) Local element index.
[in,out]connectivity_length(int *) On input, maximum length of connectivity. On output, actual length.
[out]element_connectivity(int*) The connectivity array to store connectivity in MOAB canonical numbering scheme. Array contains vertex indices in the local numbering order for vertices.
Returns
ErrCode The error code indicating success or failure.

Definition at line 1214 of file iMOAB.cpp.

1218 { 1219  assert( elem_index ); // ensure element index argument is not null 1220  assert( connectivity_length ); // ensure connectivity length argument is not null 1221  1222  appData& data = context.appDatas[*pid]; 1223  assert( ( *elem_index >= 0 ) && ( *elem_index < (int)data.primary_elems.size() ) ); 1224  1225  int num_nodes; 1226  const EntityHandle* conn; 1227  1228  EntityHandle eh = data.primary_elems[*elem_index]; 1229  1230  ErrorCode rval = context.MBI->get_connectivity( eh, conn, num_nodes );MB_CHK_ERR( rval ); 1231  1232  if( *connectivity_length < num_nodes ) 1233  { 1234  return moab::MB_FAILURE; 1235  } // wrong number of vertices 1236  1237  for( int i = 0; i < num_nodes; i++ ) 1238  { 1239  int index = data.all_verts.index( conn[i] ); 1240  1241  if( -1 == index ) 1242  { 1243  return moab::MB_FAILURE; 1244  } 1245  1246  element_connectivity[i] = index; 1247  } 1248  1249  *connectivity_length = num_nodes; 1250  1251  return moab::MB_SUCCESS; 1252 }

References appData::all_verts, GlobalContext::appDatas, context, ErrorCode, moab::Interface::get_connectivity(), moab::Range::index(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, appData::primary_elems, and moab::Range::size().

◆ iMOAB_GetElementID()

ErrCode iMOAB_GetElementID ( iMOAB_AppID  pid,
iMOAB_GlobalID global_block_ID,
int *  num_elements_in_block,
iMOAB_GlobalID global_element_ID,
iMOAB_LocalID local_element_ID 
)

Get the global IDs for all locally visible elements belonging to a particular block.

Note
The method will return also the local index of each element, in the local range that contains all visible elements.

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]global_block_ID(iMOAB_GlobalID*) The global block ID of the set being queried.
[in]num_elements_in_block(int*) The allocated size of global element ID array, same as num_elements_in_block returned from GetBlockInfo().
[out]global_element_ID(iMOAB_GlobalID*) The global IDs for all locally visible elements (array allocated by user).
[out]local_element_ID(iMOAB_LocalID*) (Optional) The local IDs for all locally visible elements (index in the range of all primary elements in the rank)
Returns
ErrCode The error code indicating success or failure.

Definition at line 1304 of file iMOAB.cpp.

1309 { 1310  assert( global_block_ID ); // ensure global block ID argument is not null 1311  assert( num_elements_in_block ); // ensure number of elements in block argument is not null 1312  1313  appData& data = context.appDatas[*pid]; 1314  std::map< int, int >& matMap = data.matIndex; 1315  1316  std::map< int, int >::iterator it = matMap.find( *global_block_ID ); 1317  1318  if( it == matMap.end() ) 1319  { 1320  return moab::MB_FAILURE; 1321  } // error in finding block with id 1322  1323  int blockIndex = matMap[*global_block_ID]; 1324  EntityHandle matMeshSet = data.mat_sets[blockIndex]; 1325  Range elems; 1326  ErrorCode rval = context.MBI->get_entities_by_handle( matMeshSet, elems );MB_CHK_ERR( rval ); 1327  1328  if( elems.empty() ) 1329  { 1330  return moab::MB_FAILURE; 1331  } 1332  1333  if( *num_elements_in_block != (int)elems.size() ) 1334  { 1335  return moab::MB_FAILURE; 1336  } // bad memory allocation 1337  1338  rval = context.MBI->tag_get_data( context.globalID_tag, elems, global_element_ID );MB_CHK_ERR( rval ); 1339  1340  // check that elems are among primary_elems in data 1341  for( int i = 0; i < *num_elements_in_block; i++ ) 1342  { 1343  local_element_ID[i] = data.primary_elems.index( elems[i] ); 1344  1345  if( -1 == local_element_ID[i] ) 1346  { 1347  return moab::MB_FAILURE; 1348  } // error, not in local primary elements 1349  } 1350  1351  return moab::MB_SUCCESS; 1352 }

References GlobalContext::appDatas, context, moab::Range::empty(), ErrorCode, moab::Interface::get_entities_by_handle(), GlobalContext::globalID_tag, moab::Range::index(), appData::mat_sets, appData::matIndex, MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, appData::primary_elems, moab::Range::size(), and moab::Interface::tag_get_data().

◆ iMOAB_GetElementOwnership()

ErrCode iMOAB_GetElementOwnership ( iMOAB_AppID  pid,
iMOAB_GlobalID global_block_ID,
int *  num_elements_in_block,
int *  element_ownership 
)

Get the element ownership within a certain block i.e., processor ID of the element owner.

Note
: Should we use local block index for input, instead of the global block ID ?

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]global_block_ID(iMOAB_GlobalID) The global block ID of the set being queried.
[in]num_elements_in_block(int*) The allocated size of ownership array, same as num_elements_in_block returned from GetBlockInfo().
[out]element_ownership(int*) The ownership array to store processor ID for all elements (array allocated by user).
Returns
ErrCode The error code indicating success or failure.

Definition at line 1254 of file iMOAB.cpp.

1258 { 1259  assert( global_block_ID ); // ensure global block ID argument is not null 1260  assert( num_elements_in_block ); // ensure number of elements in block argument is not null 1261  1262  std::map< int, int >& matMap = context.appDatas[*pid].matIndex; 1263  1264  std::map< int, int >::iterator it = matMap.find( *global_block_ID ); 1265  1266  if( it == matMap.end() ) 1267  { 1268  return moab::MB_FAILURE; 1269  } // error in finding block with id 1270  1271  int blockIndex = matMap[*global_block_ID]; 1272  EntityHandle matMeshSet = context.appDatas[*pid].mat_sets[blockIndex]; 1273  Range elems; 1274  1275  ErrorCode rval = context.MBI->get_entities_by_handle( matMeshSet, elems );MB_CHK_ERR( rval ); 1276  1277  if( elems.empty() ) 1278  { 1279  return moab::MB_FAILURE; 1280  } 1281  1282  if( *num_elements_in_block != (int)elems.size() ) 1283  { 1284  return moab::MB_FAILURE; 1285  } // bad memory allocation 1286  1287  int i = 0; 1288 #ifdef MOAB_HAVE_MPI 1289  ParallelComm* pco = context.pcomms[*pid]; 1290 #endif 1291  1292  for( Range::iterator vit = elems.begin(); vit != elems.end(); vit++, i++ ) 1293  { 1294 #ifdef MOAB_HAVE_MPI 1295  rval = pco->get_owner( *vit, element_ownership[i] );MB_CHK_ERR( rval ); 1296 #else 1297  element_ownership[i] = 0; /* owned by 0 */ 1298 #endif 1299  } 1300  1301  return moab::MB_SUCCESS; 1302 }

References GlobalContext::appDatas, moab::Range::begin(), context, moab::Range::empty(), moab::Range::end(), ErrorCode, moab::Interface::get_entities_by_handle(), moab::ParallelComm::get_owner(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, GlobalContext::pcomms, and moab::Range::size().

◆ iMOAB_GetGlobalInfo()

ErrCode iMOAB_GetGlobalInfo ( iMOAB_AppID  pid,
int *  num_global_verts,
int *  num_global_elems 
)

Get global information about number of vertices and number of elements.

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]num_global_verts(int*) The number of total vertices in the mesh.
[in]num_global_elems(MPI_Comm) The number of total elements in the mesh.
Returns
ErrCode The error code indicating success or failure.

Definition at line 2316 of file iMOAB.cpp.

2317 { 2318  appData& data = context.appDatas[*pid]; 2319  if( NULL != num_global_verts ) 2320  { 2321  *num_global_verts = data.num_global_vertices; 2322  } 2323  if( NULL != num_global_elems ) 2324  { 2325  *num_global_elems = data.num_global_elements; 2326  } 2327  2328  return moab::MB_SUCCESS; 2329 }

References GlobalContext::appDatas, context, MB_SUCCESS, appData::num_global_elements, and appData::num_global_vertices.

◆ iMOAB_GetIntTagStorage()

ErrCode iMOAB_GetIntTagStorage ( iMOAB_AppID  pid,
const iMOAB_String  tag_storage_name,
int *  num_tag_storage_length,
int *  entity_type,
int *  tag_storage_data 
)

Get the specified values in a MOAB integer Tag.

Note
Operations: Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]tag_storage_name(iMOAB_String) The tag name to store/retreive the data in MOAB.
[in]num_tag_storage_length(int*) The size of tag storage data (e.g., num_visible_vertices*components_per_entity or num_visible_elements*components_per_entity).
[in]entity_type(int*) Type=0 for vertices, and Type=1 for primary elements.
[out]tag_storage_data(int*) The array data of type int to be copied from the internal tag memory; The data is assumed to be contiguous over the local set of visible entities (either vertices or elements).
Returns
ErrCode The error code indicating success or failure.

Definition at line 1663 of file iMOAB.cpp.

1668 { 1669  ErrorCode rval; 1670  std::string tag_name( tag_storage_name ); 1671  1672  appData& data = context.appDatas[*pid]; 1673  1674  if( data.tagMap.find( tag_name ) == data.tagMap.end() ) 1675  { 1676  return moab::MB_FAILURE; 1677  } // tag not defined 1678  1679  Tag tag = data.tagMap[tag_name]; 1680  1681  int tagLength = 0; 1682  rval = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval ); 1683  1684  DataType dtype; 1685  rval = context.MBI->tag_get_data_type( tag, dtype );MB_CHK_ERR( rval ); 1686  1687  if( dtype != MB_TYPE_INTEGER ) 1688  { 1689  return moab::MB_FAILURE; 1690  } 1691  1692  // set it on a subset of entities, based on type and length 1693  Range* ents_to_get; 1694  1695  if( *ent_type == 0 ) // vertices 1696  { 1697  ents_to_get = &data.all_verts; 1698  } 1699  else // if (*ent_type == 1) 1700  { 1701  ents_to_get = &data.primary_elems; 1702  } 1703  1704  int nents_to_get = *num_tag_storage_length / tagLength; 1705  1706  if( nents_to_get > (int)ents_to_get->size() ) 1707  { 1708  return moab::MB_FAILURE; 1709  } // to many entities to get, or too little 1710  1711  // restrict the range; everything is contiguous; or not? 1712  // Range contig_range ( * ( ents_to_get->begin() ), * ( ents_to_get->begin() + nents_to_get - 1 1713  // ) ); 1714  1715  rval = context.MBI->tag_get_data( tag, *ents_to_get, tag_storage_data );MB_CHK_ERR( rval ); 1716  1717  return moab::MB_SUCCESS; // no error 1718 }

References appData::all_verts, GlobalContext::appDatas, context, ErrorCode, MB_CHK_ERR, MB_SUCCESS, MB_TYPE_INTEGER, GlobalContext::MBI, appData::primary_elems, moab::Range::size(), moab::Interface::tag_get_data(), moab::Interface::tag_get_data_type(), moab::Interface::tag_get_length(), and appData::tagMap.

◆ iMOAB_GetMeshInfo()

ErrCode iMOAB_GetMeshInfo ( iMOAB_AppID  pid,
int *  num_visible_vertices,
int *  num_visible_elements,
int *  num_visible_blocks,
int *  num_visible_surfaceBC,
int *  num_visible_vertexBC 
)

Retrieve all the important mesh information related to vertices, elements, vertex and surface boundary conditions.

Note
Number of visible vertices and cells include ghost entities. All arrays returned have size 3. Local entities are first, then ghost entities are next. Shared vertices can be owned in MOAB sense by different tasks. Ghost vertices and cells are always owned by other tasks.

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[out]num_visible_vertices(int*) The number of vertices in the current partition/process arranged. as: owned/shared only, ghosted, total_visible (array allocated by user, size := 3).
[out]num_visible_elements(int*) The number of elements in current partition/process arranged as: owned only, ghosted/shared, total_visible (array allocated by user, size := 3).
[out]num_visible_blocks(int*) The number of material sets in local mesh in current partition/process arranged as: owned only, ghosted/shared, total_visible (array allocated by user, size := 3).
[out]num_visible_surfaceBC(int*) The number of mesh surfaces that have a NEUMANN_SET BC defined in local mesh in current partition/process arranged as: owned only, ghosted/shared, total_visible (array allocated by user, size := 3).
[out]num_visible_vertexBC(int*) The number of vertices that have a DIRICHLET_SET BC defined in local mesh in current partition/process arranged as: owned only, ghosted/shared, total_visible (array allocated by user, size := 3).
Returns
ErrCode The error code indicating success or failure.

Definition at line 882 of file iMOAB.cpp.

888 { 889  ErrorCode rval; 890  appData& data = context.appDatas[*pid]; 891  EntityHandle fileSet = data.file_set; 892  893  // this will include ghost elements 894  // first clear all data ranges; this can be called after ghosting 895  if( num_visible_elements ) 896  { 897  num_visible_elements[2] = static_cast< int >( data.primary_elems.size() ); 898  // separate ghost and local/owned primary elements 899  num_visible_elements[0] = static_cast< int >( data.owned_elems.size() ); 900  num_visible_elements[1] = static_cast< int >( data.ghost_elems.size() ); 901  } 902  if( num_visible_vertices ) 903  { 904  num_visible_vertices[2] = static_cast< int >( data.all_verts.size() ); 905  num_visible_vertices[1] = static_cast< int >( data.ghost_vertices.size() ); 906  // local are those that are not ghosts; they may include shared, not owned vertices 907  num_visible_vertices[0] = num_visible_vertices[2] - num_visible_vertices[1]; 908  } 909  910  if( num_visible_blocks ) 911  { 912  rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.material_tag ), 0, 1, 913  data.mat_sets, Interface::UNION );MB_CHK_ERR( rval ); 914  915  num_visible_blocks[2] = data.mat_sets.size(); 916  num_visible_blocks[0] = num_visible_blocks[2]; 917  num_visible_blocks[1] = 0; 918  } 919  920  if( num_visible_surfaceBC ) 921  { 922  rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.neumann_tag ), 0, 1, 923  data.neu_sets, Interface::UNION );MB_CHK_ERR( rval ); 924  925  num_visible_surfaceBC[2] = 0; 926  // count how many faces are in each neu set, and how many regions are 927  // adjacent to them; 928  int numNeuSets = (int)data.neu_sets.size(); 929  930  for( int i = 0; i < numNeuSets; i++ ) 931  { 932  Range subents; 933  EntityHandle nset = data.neu_sets[i]; 934  rval = context.MBI->get_entities_by_dimension( nset, data.dimension - 1, subents );MB_CHK_ERR( rval ); 935  936  for( Range::iterator it = subents.begin(); it != subents.end(); ++it ) 937  { 938  EntityHandle subent = *it; 939  Range adjPrimaryEnts; 940  rval = context.MBI->get_adjacencies( &subent, 1, data.dimension, false, adjPrimaryEnts );MB_CHK_ERR( rval ); 941  942  num_visible_surfaceBC[2] += (int)adjPrimaryEnts.size(); 943  } 944  } 945  946  num_visible_surfaceBC[0] = num_visible_surfaceBC[2]; 947  num_visible_surfaceBC[1] = 0; 948  } 949  950  if( num_visible_vertexBC ) 951  { 952  rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.dirichlet_tag ), 0, 1, 953  data.diri_sets, Interface::UNION );MB_CHK_ERR( rval ); 954  955  num_visible_vertexBC[2] = 0; 956  int numDiriSets = (int)data.diri_sets.size(); 957  958  for( int i = 0; i < numDiriSets; i++ ) 959  { 960  Range verts; 961  EntityHandle diset = data.diri_sets[i]; 962  rval = context.MBI->get_entities_by_dimension( diset, 0, verts );MB_CHK_ERR( rval ); 963  964  num_visible_vertexBC[2] += (int)verts.size(); 965  } 966  967  num_visible_vertexBC[0] = num_visible_vertexBC[2]; 968  num_visible_vertexBC[1] = 0; 969  } 970  971  return moab::MB_SUCCESS; 972 }

References appData::all_verts, GlobalContext::appDatas, moab::Range::begin(), context, appData::dimension, appData::diri_sets, GlobalContext::dirichlet_tag, moab::Range::end(), ErrorCode, appData::file_set, moab::Interface::get_adjacencies(), moab::Interface::get_entities_by_dimension(), moab::Interface::get_entities_by_type_and_tag(), appData::ghost_elems, appData::ghost_vertices, appData::mat_sets, GlobalContext::material_tag, MB_CHK_ERR, MB_SUCCESS, MBENTITYSET, GlobalContext::MBI, appData::neu_sets, GlobalContext::neumann_tag, appData::owned_elems, appData::primary_elems, moab::Range::size(), and moab::Interface::UNION.

◆ iMOAB_GetNeighborElements()

ErrCode iMOAB_GetNeighborElements ( iMOAB_AppID  pid,
iMOAB_LocalID local_index,
int *  num_adjacent_elements,
iMOAB_LocalID adjacent_element_IDs 
)

Retrieve the adjacencies for the element entities.

Note
Operations: Not Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]local_index(iMOAB_LocalID*) The local element ID for which adjacency information is needed.
[out]num_adjacent_elements(int*) The total number of adjacent elements.
[out]adjacent_element_IDs(iMOAB_LocalID*) The local element IDs of all adjacent elements to the current one (typically, num_total_sides for internal elements or num_total_sides-num_sides_on_boundary for boundary elements).
Returns
ErrCode The error code indicating success or failure.

Definition at line 2184 of file iMOAB.cpp.

2188 { 2189  ErrorCode rval; 2190  2191  // one neighbor for each subentity of dimension-1 2192  MeshTopoUtil mtu( context.MBI ); 2193  appData& data = context.appDatas[*pid]; 2194  EntityHandle eh = data.primary_elems[*local_index]; 2195  Range adjs; 2196  rval = mtu.get_bridge_adjacencies( eh, data.dimension - 1, data.dimension, adjs );MB_CHK_ERR( rval ); 2197  2198  if( *num_adjacent_elements < (int)adjs.size() ) 2199  { 2200  return moab::MB_FAILURE; 2201  } // not dimensioned correctly 2202  2203  *num_adjacent_elements = (int)adjs.size(); 2204  2205  for( int i = 0; i < *num_adjacent_elements; i++ ) 2206  { 2207  adjacent_element_IDs[i] = data.primary_elems.index( adjs[i] ); 2208  } 2209  2210  return moab::MB_SUCCESS; 2211 }

References GlobalContext::appDatas, context, appData::dimension, ErrorCode, moab::MeshTopoUtil::get_bridge_adjacencies(), moab::Range::index(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, appData::primary_elems, and moab::Range::size().

◆ iMOAB_GetNeighborVertices()

ErrCode iMOAB_GetNeighborVertices ( iMOAB_AppID  pid,
iMOAB_LocalID local_vertex_ID,
int *  num_adjacent_vertices,
iMOAB_LocalID adjacent_vertex_IDs 
)

Get the adjacencies for the vertex entities.

Note
Operations: Not Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]local_vertex_ID(iMOAB_LocalID*) The local vertex ID for which adjacency information is needed.
[out]num_adjacent_vertices(int*) The total number of adjacent vertices.
[out]adjacent_vertex_IDs(iMOAB_LocalID*) The local element IDs of all adjacent vertices to the current one (typically, num_total_sides for internal elements or num_total_sides-num_sides_on_boundary for boundary elements).
Returns
ErrCode The error code indicating success or failure.

◆ iMOAB_GetPointerToSurfaceBC()

ErrCode iMOAB_GetPointerToSurfaceBC ( iMOAB_AppID  pid,
int *  surface_BC_length,
iMOAB_LocalID local_element_ID,
int *  reference_surface_ID,
int *  boundary_condition_value 
)

Get the surface boundary condition information.

Note
Operations: Not Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]surface_BC_length(int*) The allocated size of surface boundary condition array, same as num_visible_surfaceBC returned by GetMeshInfo().
[out]local_element_ID(iMOAB_LocalID*) The local element IDs that contains the side with the surface BC.
[out]reference_surface_ID(int*) The surface number with the BC in the canonical reference element (e.g., 1 to 6 for HEX, 1-4 for TET).
[out]boundary_condition_value(int*) The boundary condition type as obtained from the mesh description (value of the NeumannSet defined on the element).
Returns
ErrCode The error code indicating success or failure.

Definition at line 1354 of file iMOAB.cpp.

1359 { 1360  // we have to fill bc data for neumann sets;/ 1361  ErrorCode rval; 1362  1363  // it was counted above, in GetMeshInfo 1364  appData& data = context.appDatas[*pid]; 1365  int numNeuSets = (int)data.neu_sets.size(); 1366  1367  int index = 0; // index [0, surface_BC_length) for the arrays returned 1368  1369  for( int i = 0; i < numNeuSets; i++ ) 1370  { 1371  Range subents; 1372  EntityHandle nset = data.neu_sets[i]; 1373  rval = context.MBI->get_entities_by_dimension( nset, data.dimension - 1, subents );MB_CHK_ERR( rval ); 1374  1375  int neuVal; 1376  rval = context.MBI->tag_get_data( context.neumann_tag, &nset, 1, &neuVal );MB_CHK_ERR( rval ); 1377  1378  for( Range::iterator it = subents.begin(); it != subents.end(); ++it ) 1379  { 1380  EntityHandle subent = *it; 1381  Range adjPrimaryEnts; 1382  rval = context.MBI->get_adjacencies( &subent, 1, data.dimension, false, adjPrimaryEnts );MB_CHK_ERR( rval ); 1383  1384  // get global id of the primary ents, and side number of the quad/subentity 1385  // this is moab ordering 1386  for( Range::iterator pit = adjPrimaryEnts.begin(); pit != adjPrimaryEnts.end(); pit++ ) 1387  { 1388  EntityHandle primaryEnt = *pit; 1389  // get global id 1390  /*int globalID; 1391  rval = context.MBI->tag_get_data(gtags[3], &primaryEnt, 1, &globalID); 1392  if (MB_SUCCESS!=rval) 1393  return moab::MB_FAILURE; 1394  global_element_ID[index] = globalID;*/ 1395  1396  // get local element id 1397  local_element_ID[index] = data.primary_elems.index( primaryEnt ); 1398  1399  if( -1 == local_element_ID[index] ) 1400  { 1401  return moab::MB_FAILURE; 1402  } // did not find the element locally 1403  1404  int side_number, sense, offset; 1405  rval = context.MBI->side_number( primaryEnt, subent, side_number, sense, offset );MB_CHK_ERR( rval ); 1406  1407  reference_surface_ID[index] = side_number + 1; // moab is from 0 to 5, it needs 1 to 6 1408  boundary_condition_value[index] = neuVal; 1409  index++; 1410  } 1411  } 1412  } 1413  1414  if( index != *surface_BC_length ) 1415  { 1416  return moab::MB_FAILURE; 1417  } // error in array allocations 1418  1419  return moab::MB_SUCCESS; 1420 }

References GlobalContext::appDatas, moab::Range::begin(), context, appData::dimension, moab::Range::end(), ErrorCode, moab::Interface::get_adjacencies(), moab::Interface::get_entities_by_dimension(), moab::Range::index(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, appData::neu_sets, GlobalContext::neumann_tag, appData::primary_elems, moab::Interface::side_number(), moab::side_number(), moab::Range::size(), and moab::Interface::tag_get_data().

◆ iMOAB_GetPointerToVertexBC()

ErrCode iMOAB_GetPointerToVertexBC ( iMOAB_AppID  pid,
int *  vertex_BC_length,
iMOAB_LocalID local_vertex_ID,
int *  boundary_condition_value 
)

Get the vertex boundary condition information.

Note
Operations: Not Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]vertex_BC_length(int) The allocated size of vertex boundary condition array, same as num_visible_vertexBC returned by GetMeshInfo().
[out]local_vertex_ID(iMOAB_LocalID*) The local vertex ID that has Dirichlet BC defined.
[out]boundary_condition_value(int*) The boundary condition type as obtained from the mesh description (value of the Dirichlet_Set tag defined on the vertex).
Returns
ErrCode The error code indicating success or failure.

Definition at line 1422 of file iMOAB.cpp.

1426 { 1427  ErrorCode rval; 1428  1429  // it was counted above, in GetMeshInfo 1430  appData& data = context.appDatas[*pid]; 1431  int numDiriSets = (int)data.diri_sets.size(); 1432  int index = 0; // index [0, *vertex_BC_length) for the arrays returned 1433  1434  for( int i = 0; i < numDiriSets; i++ ) 1435  { 1436  Range verts; 1437  EntityHandle diset = data.diri_sets[i]; 1438  rval = context.MBI->get_entities_by_dimension( diset, 0, verts );MB_CHK_ERR( rval ); 1439  1440  int diriVal; 1441  rval = context.MBI->tag_get_data( context.dirichlet_tag, &diset, 1, &diriVal );MB_CHK_ERR( rval ); 1442  1443  for( Range::iterator vit = verts.begin(); vit != verts.end(); ++vit ) 1444  { 1445  EntityHandle vt = *vit; 1446  /*int vgid; 1447  rval = context.MBI->tag_get_data(gtags[3], &vt, 1, &vgid); 1448  if (MB_SUCCESS!=rval) 1449  return moab::MB_FAILURE; 1450  global_vertext_ID[index] = vgid;*/ 1451  local_vertex_ID[index] = data.all_verts.index( vt ); 1452  1453  if( -1 == local_vertex_ID[index] ) 1454  { 1455  return moab::MB_FAILURE; 1456  } // vertex was not found 1457  1458  boundary_condition_value[index] = diriVal; 1459  index++; 1460  } 1461  } 1462  1463  if( *vertex_BC_length != index ) 1464  { 1465  return moab::MB_FAILURE; 1466  } // array allocation issue 1467  1468  return moab::MB_SUCCESS; 1469 }

References appData::all_verts, GlobalContext::appDatas, moab::Range::begin(), context, appData::diri_sets, GlobalContext::dirichlet_tag, moab::Range::end(), ErrorCode, moab::Interface::get_entities_by_dimension(), moab::Range::index(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, moab::Range::size(), and moab::Interface::tag_get_data().

◆ iMOAB_GetVertexID()

ErrCode iMOAB_GetVertexID ( iMOAB_AppID  pid,
int *  vertices_length,
iMOAB_GlobalID global_vertex_ID 
)

Get the global vertex IDs for all locally visible (owned and shared/ghosted) vertices.

Note
The array should be allocated by the user, sized with the total number of visible vertices from iMOAB_GetMeshInfo method.

Operations: Not collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]vertices_length(int*) The allocated size of array (typical size := num_visible_vertices).
[out]global_vertex_ID(iMOAB_GlobalID*) The global IDs for all locally visible vertices (array allocated by user).
Returns
ErrCode The error code indicating success or failure.

Definition at line 974 of file iMOAB.cpp.

975 { 976  IMOAB_CHECKPOINTER( vertices_length, 2 ); 977  IMOAB_CHECKPOINTER( global_vertex_ID, 3 ); 978  979  const Range& verts = context.appDatas[*pid].all_verts; 980  // check for problems with array length 981  IMOAB_ASSERT( *vertices_length == static_cast< int >( verts.size() ), "Invalid vertices length provided" ); 982  983  // global id tag is context.globalID_tag 984  return context.MBI->tag_get_data( context.globalID_tag, verts, global_vertex_ID ); 985 }

References GlobalContext::appDatas, context, GlobalContext::globalID_tag, IMOAB_ASSERT, IMOAB_CHECKPOINTER, GlobalContext::MBI, moab::Range::size(), and moab::Interface::tag_get_data().

◆ iMOAB_GetVertexOwnership()

ErrCode iMOAB_GetVertexOwnership ( iMOAB_AppID  pid,
int *  vertices_length,
int *  visible_global_rank_ID 
)

Get vertex ownership information.

For each vertex based on the local ID, return the process that owns the vertex (local, shared or ghost)

Note
Shared vertices could be owned by different tasks. Local and shared vertices are first, ghost vertices are next in the array. Ghost vertices are always owned by a different process ID. Array allocated by the user with total size of visible vertices.

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]vertices_length(int*) The allocated size of array (typically size := num_visible_vertices).
[out]visible_global_rank_ID(int*) The processor rank owning each of the local vertices.
Returns
ErrCode The error code indicating success or failure.

Definition at line 987 of file iMOAB.cpp.

988 { 989  assert( vertices_length && *vertices_length ); 990  assert( visible_global_rank_ID ); 991  992  Range& verts = context.appDatas[*pid].all_verts; 993  int i = 0; 994 #ifdef MOAB_HAVE_MPI 995  ParallelComm* pco = context.pcomms[*pid]; 996  997  for( Range::iterator vit = verts.begin(); vit != verts.end(); vit++, i++ ) 998  { 999  ErrorCode rval = pco->get_owner( *vit, visible_global_rank_ID[i] );MB_CHK_ERR( rval ); 1000  } 1001  1002  if( i != *vertices_length ) 1003  { 1004  return moab::MB_FAILURE; 1005  } // warning array allocation problem 1006  1007 #else 1008  1009  /* everything owned by proc 0 */ 1010  if( (int)verts.size() != *vertices_length ) 1011  { 1012  return moab::MB_FAILURE; 1013  } // warning array allocation problem 1014  1015  for( Range::iterator vit = verts.begin(); vit != verts.end(); vit++, i++ ) 1016  { 1017  visible_global_rank_ID[i] = 0; 1018  } // all vertices are owned by processor 0, as this is serial run 1019  1020 #endif 1021  1022  return moab::MB_SUCCESS; 1023 }

References GlobalContext::appDatas, moab::Range::begin(), context, moab::Range::end(), ErrorCode, moab::ParallelComm::get_owner(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::pcomms, and moab::Range::size().

◆ iMOAB_GetVisibleElementsInfo()

ErrCode iMOAB_GetVisibleElementsInfo ( iMOAB_AppID  pid,
int *  num_visible_elements,
iMOAB_GlobalID element_global_IDs,
int *  ranks,
iMOAB_GlobalID block_IDs 
)

Get the visible elements information.

Return for all visible elements the global IDs, ranks they belong to, block ids they belong to.

Note
Operations: Not Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]num_visible_elements(int*) The number of visible elements (returned by GetMeshInfo).
[out]element_global_IDs(iMOAB_GlobalID*) Element global ids to be added to the block.
[out]ranks(int*) The owning ranks of elements.
[out]block_IDs(iMOAB_GlobalID*) The block ids containing the elements.
Returns
ErrCode The error code indicating success or failure.

Definition at line 1106 of file iMOAB.cpp.

1111 { 1112  appData& data = context.appDatas[*pid]; 1113 #ifdef MOAB_HAVE_MPI 1114  ParallelComm* pco = context.pcomms[*pid]; 1115 #endif 1116  1117  ErrorCode rval = context.MBI->tag_get_data( context.globalID_tag, data.primary_elems, element_global_IDs );MB_CHK_ERR( rval ); 1118  1119  int i = 0; 1120  1121  for( Range::iterator eit = data.primary_elems.begin(); eit != data.primary_elems.end(); ++eit, ++i ) 1122  { 1123 #ifdef MOAB_HAVE_MPI 1124  rval = pco->get_owner( *eit, ranks[i] );MB_CHK_ERR( rval ); 1125  1126 #else 1127  /* everything owned by task 0 */ 1128  ranks[i] = 0; 1129 #endif 1130  } 1131  1132  for( Range::iterator mit = data.mat_sets.begin(); mit != data.mat_sets.end(); ++mit ) 1133  { 1134  EntityHandle matMeshSet = *mit; 1135  Range elems; 1136  rval = context.MBI->get_entities_by_handle( matMeshSet, elems );MB_CHK_ERR( rval ); 1137  1138  int valMatTag; 1139  rval = context.MBI->tag_get_data( context.material_tag, &matMeshSet, 1, &valMatTag );MB_CHK_ERR( rval ); 1140  1141  for( Range::iterator eit = elems.begin(); eit != elems.end(); ++eit ) 1142  { 1143  EntityHandle eh = *eit; 1144  int index = data.primary_elems.index( eh ); 1145  1146  if( -1 == index ) 1147  { 1148  return moab::MB_FAILURE; 1149  } 1150  1151  if( -1 >= *num_visible_elements ) 1152  { 1153  return moab::MB_FAILURE; 1154  } 1155  1156  block_IDs[index] = valMatTag; 1157  } 1158  } 1159  1160  return moab::MB_SUCCESS; 1161 }

References GlobalContext::appDatas, moab::Range::begin(), context, moab::Range::end(), ErrorCode, moab::Interface::get_entities_by_handle(), moab::ParallelComm::get_owner(), GlobalContext::globalID_tag, moab::Range::index(), appData::mat_sets, GlobalContext::material_tag, MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, GlobalContext::pcomms, appData::primary_elems, and moab::Interface::tag_get_data().

◆ iMOAB_GetVisibleVerticesCoordinates()

ErrCode iMOAB_GetVisibleVerticesCoordinates ( iMOAB_AppID  pid,
int *  coords_length,
double *  coordinates 
)

Get vertex coordinates for all local (owned and ghosted) vertices.

Note
coordinates are returned in an array allocated by user, interleaved. (do need an option for blocked coordinates ?) size of the array is dimension times number of visible vertices. The local ordering is implicit, owned/shared vertices are first, then ghosts.

Operations: Not Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]coords_length(int*) The size of the allocated coordinate array (array allocated by user, size := 3*num_visible_vertices).
[out]coordinates(double*) The pointer to user allocated memory that will be filled with interleaved coordinates.
Returns
ErrCode The error code indicating success or failure.

Definition at line 1025 of file iMOAB.cpp.

1026 { 1027  Range& verts = context.appDatas[*pid].all_verts; 1028  1029  // interleaved coordinates, so that means deep copy anyway 1030  if( *coords_length != 3 * (int)verts.size() ) 1031  { 1032  return moab::MB_FAILURE; 1033  } 1034  1035  ErrorCode rval = context.MBI->get_coords( verts, coordinates );MB_CHK_ERR( rval ); 1036  1037  return moab::MB_SUCCESS; 1038 }

References GlobalContext::appDatas, context, ErrorCode, moab::Interface::get_coords(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, and moab::Range::size().

◆ iMOAB_Initialize()

ErrCode iMOAB_Initialize ( int  argc,
iMOAB_String argv 
)

Initialize the iMOAB interface implementation.

Note
Will create the MOAB instance, if not created already (reference counted).

Operations: Collective

Parameters
[in]argc(int) Number of command line arguments
[in]argv(iMOAB_String*) Command line arguments
Returns
ErrCode

Definition at line 138 of file iMOAB.cpp.

139 { 140  if( argc ) IMOAB_CHECKPOINTER( argv, 1 ); 141  142  context.iArgc = argc; 143  context.iArgv = argv; // shallow copy 144  145  if( 0 == context.refCountMB ) 146  { 147  context.MBI = new( std::nothrow ) moab::Core; 148  // retrieve the default tags 149  const char* const shared_set_tag_names[] = { MATERIAL_SET_TAG_NAME, NEUMANN_SET_TAG_NAME, 150  DIRICHLET_SET_TAG_NAME, GLOBAL_ID_TAG_NAME }; 151  // blocks, visible surfaceBC(neumann), vertexBC (Dirichlet), global id, parallel partition 152  Tag gtags[4]; 153  154  for( int i = 0; i < 4; i++ ) 155  { 156  157  ErrorCode rval = 158  context.MBI->tag_get_handle( shared_set_tag_names[i], 1, MB_TYPE_INTEGER, gtags[i], MB_TAG_ANY );MB_CHK_ERR( rval ); 159  } 160  161  context.material_tag = gtags[0]; 162  context.neumann_tag = gtags[1]; 163  context.dirichlet_tag = gtags[2]; 164  context.globalID_tag = gtags[3]; 165  } 166  167  context.MPI_initialized = false; 168 #ifdef MOAB_HAVE_MPI 169  int flagInit; 170  MPI_Initialized( &flagInit ); 171  172  if( flagInit && !context.MPI_initialized ) 173  { 174  MPI_Comm_size( MPI_COMM_WORLD, &context.worldprocs ); 175  MPI_Comm_rank( MPI_COMM_WORLD, &context.globalrank ); 176  context.MPI_initialized = true; 177  } 178 #endif 179  180  context.refCountMB++; 181  return moab::MB_SUCCESS; 182 }

References context, DIRICHLET_SET_TAG_NAME, GlobalContext::dirichlet_tag, ErrorCode, GLOBAL_ID_TAG_NAME, GlobalContext::globalID_tag, GlobalContext::globalrank, GlobalContext::iArgc, GlobalContext::iArgv, IMOAB_CHECKPOINTER, MATERIAL_SET_TAG_NAME, GlobalContext::material_tag, MB_CHK_ERR, MB_SUCCESS, MB_TAG_ANY, MB_TYPE_INTEGER, GlobalContext::MBI, GlobalContext::MPI_initialized, NEUMANN_SET_TAG_NAME, GlobalContext::neumann_tag, GlobalContext::refCountMB, moab::Interface::tag_get_handle(), and GlobalContext::worldprocs.

Referenced by iMOAB_InitializeFortran().

◆ iMOAB_InitializeFortran()

ErrCode iMOAB_InitializeFortran ( void  )

Initialize the iMOAB interface implementation from Fortran driver.

It will create the MOAB instance, if not created already (reference counted).

Operations: Collective

Returns
ErrCode The error code indicating success or failure.

Definition at line 184 of file iMOAB.cpp.

185 { 186  return iMOAB_Initialize( 0, 0 ); 187 }

References iMOAB_Initialize().

◆ iMOAB_LoadMesh()

ErrCode iMOAB_LoadMesh ( iMOAB_AppID  pid,
const iMOAB_String  filename,
const iMOAB_String  read_options,
int *  num_ghost_layers 
)

Load a MOAB mesh file in parallel and exchange ghost layers as requested.

Note
All communication is MPI-based, and read options include parallel loading information, resolving shared entities. Local MOAB instance is populated with mesh cells and vertices in the corresponding local partitions.

This will also exchange ghost cells and vertices, as requested. The default bridge dimension is 0 (vertices), and all additional lower dimensional sub-entities are exchanged (mesh edges and faces). The tags in the file are not exchanged by default. Default tag information for GLOBAL_ID, MATERIAL_SET, NEUMANN_SET and DIRICHLET_SET is exchanged. Global ID tag is exchanged for all cells and vertices. Material sets, Neumann sets and Dirichlet sets are all augmented with the ghost entities.

Operations: Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]filename(iMOAB_String) The MOAB mesh file (H5M) to load onto the internal application mesh set.
[in]read_options(iMOAB_String) Additional options for reading the MOAB mesh file in parallel.
[in]num_ghost_layers(int*) The total number of ghost layers to exchange during mesh loading.
Returns
ErrCode The error code indicating success or failure.

Definition at line 611 of file iMOAB.cpp.

615 { 616  IMOAB_CHECKPOINTER( filename, 2 ); 617  IMOAB_ASSERT( strlen( filename ), "Invalid filename length." ); 618  IMOAB_CHECKPOINTER( num_ghost_layers, 4 ); 619  620  // make sure we use the file set and pcomm associated with the *pid 621  std::ostringstream newopts; 622  if( read_options ) newopts << read_options; 623  624 #ifdef MOAB_HAVE_MPI 625  626  if( context.MPI_initialized ) 627  { 628  if( context.worldprocs > 1 ) 629  { 630  std::string opts( ( read_options ? read_options : "" ) ); 631  std::string pcid( "PARALLEL_COMM=" ); 632  std::size_t found = opts.find( pcid ); 633  634  if( found != std::string::npos ) 635  { 636  std::cerr << " cannot specify PARALLEL_COMM option, it is implicit \n"; 637  return moab::MB_FAILURE; 638  } 639  640  // in serial, apply PARALLEL_COMM option only for h5m files; it does not work for .g 641  // files (used in test_remapping) 642  std::string filen( filename ); 643  std::string::size_type idx = filen.rfind( '.' ); 644  645  if( idx != std::string::npos ) 646  { 647  std::string extension = filen.substr( idx + 1 ); 648  if( ( extension == std::string( "h5m" ) ) || ( extension == std::string( "nc" ) ) ) 649  newopts << ";;PARALLEL_COMM=" << *pid; 650  } 651  652  if( *num_ghost_layers >= 1 ) 653  { 654  // if we want ghosts, we will want additional entities, the last .1 655  // because the addl ents can be edges, faces that are part of the neumann sets 656  std::string pcid2( "PARALLEL_GHOSTS=" ); 657  std::size_t found2 = opts.find( pcid2 ); 658  659  if( found2 != std::string::npos ) 660  { 661  std::cout << " PARALLEL_GHOSTS option is already specified, ignore passed " 662  "number of layers \n"; 663  } 664  else 665  { 666  // dimension of primary entities is 3 here, but it could be 2 for climate 667  // meshes; we would need to pass PARALLEL_GHOSTS explicitly for 2d meshes, for 668  // example: ";PARALLEL_GHOSTS=2.0.1" 669  newopts << ";PARALLEL_GHOSTS=3.0." << *num_ghost_layers << ".3"; 670  } 671  } 672  } 673  } 674 #else 675  IMOAB_ASSERT( *num_ghost_layers == 0, "Cannot provide ghost layers in serial." ); 676 #endif 677  678  // Now let us actually load the MOAB file with the appropriate read options 679  ErrorCode rval = context.MBI->load_file( filename, &context.appDatas[*pid].file_set, newopts.str().c_str() );MB_CHK_ERR( rval ); 680  681 #ifdef VERBOSE 682  // some debugging stuff 683  std::ostringstream outfile; 684 #ifdef MOAB_HAVE_MPI 685  int rank = context.pcomms[*pid]->rank(); 686  int nprocs = context.pcomms[*pid]->size(); 687  outfile << "TaskMesh_n" << nprocs << "." << rank << ".h5m"; 688 #else 689  outfile << "TaskMesh_n1.0.h5m"; 690 #endif 691  // the mesh contains ghosts too, but they are not part of mat/neumann set 692  // write in serial the file, to see what tags are missing 693  rval = context.MBI->write_file( outfile.str().c_str() );MB_CHK_ERR( rval ); // everything on current task, written in serial 694 #endif 695  696  // Update ghost layer information 697  context.appDatas[*pid].num_ghost_layers = *num_ghost_layers; 698  699  // Update mesh information 700  return iMOAB_UpdateMeshInfo( pid ); 701 }

References GlobalContext::appDatas, context, ErrorCode, IMOAB_ASSERT, IMOAB_CHECKPOINTER, iMOAB_UpdateMeshInfo(), moab::Interface::load_file(), MB_CHK_ERR, GlobalContext::MBI, GlobalContext::MPI_initialized, GlobalContext::pcomms, GlobalContext::worldprocs, and moab::Interface::write_file().

◆ iMOAB_ReadHeaderInfo()

ErrCode iMOAB_ReadHeaderInfo ( const iMOAB_String  filename,
int *  num_global_vertices,
int *  num_global_elements,
int *  num_dimension,
int *  num_parts 
)

It should be called on master task only, and information obtained could be broadcasted by the user. It is a fast lookup in the header of the file.

Note
Operations: Not collective
Parameters
[in]filename(iMOAB_String) The MOAB mesh file (H5M) to probe for header information.
[out]num_global_vertices(int*) The total number of vertices in the mesh file.
[out]num_global_elements(int*) The total number of elements (of highest dimension only).
[out]num_dimension(int*) The highest dimension of elements in the mesh (Edge=1, Tri/Quad=2, Tet/Hex/Prism/Pyramid=3).
[out]num_parts(int*) The total number of partitions available in the mesh file, typically partitioned with mbpart during pre-processing.
Returns
ErrCode The error code indicating success or failure.

Definition at line 467 of file iMOAB.cpp.

472 { 473  IMOAB_CHECKPOINTER( filename, 1 ); 474  IMOAB_ASSERT( strlen( filename ), "Invalid filename length." ); 475  476 #ifdef MOAB_HAVE_HDF5 477  std::string filen( filename ); 478  479  int edges = 0; 480  int faces = 0; 481  int regions = 0; 482  if( num_global_vertices ) *num_global_vertices = 0; 483  if( num_global_elements ) *num_global_elements = 0; 484  if( num_dimension ) *num_dimension = 0; 485  if( num_parts ) *num_parts = 0; 486  487  mhdf_FileHandle file; 488  mhdf_Status status; 489  unsigned long max_id; 490  struct mhdf_FileDesc* data; 491  492  file = mhdf_openFile( filen.c_str(), 0, &max_id, -1, &status ); 493  494  if( mhdf_isError( &status ) ) 495  { 496  fprintf( stderr, "%s: %s\n", filename, mhdf_message( &status ) ); 497  return moab::MB_FAILURE; 498  } 499  500  data = mhdf_getFileSummary( file, H5T_NATIVE_ULONG, &status, 501  1 ); // will use extra set info; will get parallel partition tag info too! 502  503  if( mhdf_isError( &status ) ) 504  { 505  fprintf( stderr, "%s: %s\n", filename, mhdf_message( &status ) ); 506  return moab::MB_FAILURE; 507  } 508  509  if( num_dimension ) *num_dimension = data->nodes.vals_per_ent; 510  if( num_global_vertices ) *num_global_vertices = (int)data->nodes.count; 511  512  for( int i = 0; i < data->num_elem_desc; i++ ) 513  { 514  struct mhdf_ElemDesc* el_desc = &( data->elems[i] ); 515  struct mhdf_EntDesc* ent_d = &( el_desc->desc ); 516  517  if( 0 == strcmp( el_desc->type, mhdf_EDGE_TYPE_NAME ) ) 518  { 519  edges += ent_d->count; 520  } 521  522  if( 0 == strcmp( el_desc->type, mhdf_TRI_TYPE_NAME ) ) 523  { 524  faces += ent_d->count; 525  } 526  527  if( 0 == strcmp( el_desc->type, mhdf_QUAD_TYPE_NAME ) ) 528  { 529  faces += ent_d->count; 530  } 531  532  if( 0 == strcmp( el_desc->type, mhdf_POLYGON_TYPE_NAME ) ) 533  { 534  faces += ent_d->count; 535  } 536  537  if( 0 == strcmp( el_desc->type, mhdf_TET_TYPE_NAME ) ) 538  { 539  regions += ent_d->count; 540  } 541  542  if( 0 == strcmp( el_desc->type, mhdf_PYRAMID_TYPE_NAME ) ) 543  { 544  regions += ent_d->count; 545  } 546  547  if( 0 == strcmp( el_desc->type, mhdf_PRISM_TYPE_NAME ) ) 548  { 549  regions += ent_d->count; 550  } 551  552  if( 0 == strcmp( el_desc->type, mdhf_KNIFE_TYPE_NAME ) ) 553  { 554  regions += ent_d->count; 555  } 556  557  if( 0 == strcmp( el_desc->type, mdhf_HEX_TYPE_NAME ) ) 558  { 559  regions += ent_d->count; 560  } 561  562  if( 0 == strcmp( el_desc->type, mhdf_POLYHEDRON_TYPE_NAME ) ) 563  { 564  regions += ent_d->count; 565  } 566  567  if( 0 == strcmp( el_desc->type, mhdf_SEPTAHEDRON_TYPE_NAME ) ) 568  { 569  regions += ent_d->count; 570  } 571  } 572  573  if( num_parts ) *num_parts = data->numEntSets[0]; 574  575  // is this required? 576  if( edges > 0 ) 577  { 578  if( num_dimension ) *num_dimension = 1; // I don't think it will ever return 1 579  if( num_global_elements ) *num_global_elements = edges; 580  } 581  582  if( faces > 0 ) 583  { 584  if( num_dimension ) *num_dimension = 2; 585  if( num_global_elements ) *num_global_elements = faces; 586  } 587  588  if( regions > 0 ) 589  { 590  if( num_dimension ) *num_dimension = 3; 591  if( num_global_elements ) *num_global_elements = regions; 592  } 593  594  mhdf_closeFile( file, &status ); 595  596  free( data ); 597  598 #else 599  std::cout << filename 600  << ": Please reconfigure with HDF5. Cannot retrieve header information for file " 601  "formats other than a h5m file.\n"; 602  if( num_global_vertices ) *num_global_vertices = 0; 603  if( num_global_elements ) *num_global_elements = 0; 604  if( num_dimension ) *num_dimension = 0; 605  if( num_parts ) *num_parts = 0; 606 #endif 607  608  return moab::MB_SUCCESS; 609 }

References mhdf_EntDesc::count, mhdf_ElemDesc::desc, mhdf_FileDesc::elems, IMOAB_ASSERT, IMOAB_CHECKPOINTER, MB_SUCCESS, mdhf_HEX_TYPE_NAME, mdhf_KNIFE_TYPE_NAME, mhdf_closeFile(), mhdf_EDGE_TYPE_NAME, mhdf_getFileSummary(), mhdf_isError(), mhdf_message(), mhdf_openFile(), mhdf_POLYGON_TYPE_NAME, mhdf_POLYHEDRON_TYPE_NAME, mhdf_PRISM_TYPE_NAME, mhdf_PYRAMID_TYPE_NAME, mhdf_QUAD_TYPE_NAME, mhdf_SEPTAHEDRON_TYPE_NAME, mhdf_TET_TYPE_NAME, mhdf_TRI_TYPE_NAME, mhdf_FileDesc::nodes, mhdf_FileDesc::num_elem_desc, mhdf_FileDesc::numEntSets, mhdf_ElemDesc::type, and mhdf_EntDesc::vals_per_ent.

◆ iMOAB_ReduceTagsMax()

ErrCode iMOAB_ReduceTagsMax ( iMOAB_AppID  pid,
int *  tag_index,
int *  ent_type 
)

Perform global reductions with the processes in the current applications communicator. Specifically this routine performs reductions on the maximum value of the tag indicated by tag_index.

Note
Operations: Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]tag_index(int*) The tag index of interest.
[in]ent_type(int*) The type of entity for tag reduction operation (vertices = 0, elements = 1)
Returns
ErrCode The error code indicating success or failure.

Definition at line 2144 of file iMOAB.cpp.

2145 { 2146 #ifdef MOAB_HAVE_MPI 2147  appData& data = context.appDatas[*pid]; 2148  Range ent_exchange; 2149  2150  if( *tag_index < 0 || *tag_index >= (int)data.tagList.size() ) 2151  { 2152  return moab::MB_FAILURE; 2153  } // error in tag index 2154  2155  Tag tagh = data.tagList[*tag_index]; 2156  2157  if( *ent_type == 0 ) 2158  { 2159  ent_exchange = data.all_verts; 2160  } 2161  else if( *ent_type == 1 ) 2162  { 2163  ent_exchange = data.primary_elems; 2164  } 2165  else 2166  { 2167  return moab::MB_FAILURE; 2168  } // unexpected type 2169  2170  ParallelComm* pco = context.pcomms[*pid]; 2171  // we could do different MPI_Op; do not bother now, we will call from fortran 2172  ErrorCode rval = pco->reduce_tags( tagh, MPI_MAX, ent_exchange );MB_CHK_ERR( rval ); 2173  2174 #else 2175  /* do nothing if serial */ 2176  // just silence the warning 2177  // do not call sync tags in serial! 2178  int k = *pid + *tag_index + *ent_type; 2179  k++; // just do junk, to avoid complaints 2180 #endif 2181  return moab::MB_SUCCESS; 2182 }

References appData::all_verts, GlobalContext::appDatas, context, ErrorCode, MB_CHK_ERR, MB_SUCCESS, GlobalContext::pcomms, appData::primary_elems, moab::ParallelComm::reduce_tags(), and appData::tagList.

◆ iMOAB_RegisterApplication()

ErrCode iMOAB_RegisterApplication ( const iMOAB_String  app_name,
int *  compid,
iMOAB_AppID  pid 
)

Register application - Create a unique application ID and bootstrap interfaces for further queries.

Note
Internally, a mesh set will be associated with the application ID and all subsequent queries on the MOAB instance will be directed to this mesh/file set.

Operations: Collective

Parameters
[in]app_name(iMOAB_String) Application name (PROTEUS, NEK5000, etc)
[in]comm(MPI_Comm*) MPI communicator to be used for all mesh-related queries originating from this application
[in]compid(int*) The unique external component identifier
[out]pid(iMOAB_AppID) The unique pointer to the application ID
Returns
ErrCode The error code indicating success or failure.

◆ iMOAB_RegisterApplicationFortran()

ErrCode iMOAB_RegisterApplicationFortran ( const iMOAB_String  app_name,
int *  compid,
iMOAB_AppID  pid 
)

Register a Fortran-based application - Create a unique application ID and bootstrap interfaces for further queries.

Note
Internally, the comm object, usually a 32 bit integer in Fortran, will be converted using MPI_Comm_f2c and stored as MPI_Comm. A mesh set will be associated with the application ID and all subsequent queries on the MOAB instance will be directed to this mesh/file set.

Operations: Collective

Parameters
[in]app_name(iMOAB_String) Application name ("MySolver", "E3SM", "DMMoab", "PROTEUS", "NEK5000", etc)
[in]communicator(int*) Fortran MPI communicator to be used for all mesh-related queries originating from this application.
[in]compid(int*) External component id, (A const unique identifier).
[out]pid(iMOAB_AppID) The unique pointer to the application ID.
Returns
ErrCode The error code indicating success or failure.

◆ iMOAB_ResolveSharedEntities()

ErrCode iMOAB_ResolveSharedEntities ( iMOAB_AppID  pid,
int *  num_verts,
int *  marker 
)

Resolve shared entities using global markers on shared vertices.

Note
Global markers can be a global node id, for example, or a global DoF number (as for HOMME)

Operations: Collective .

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]num_verts(int*) Number of vertices.
[in]marker(int*) Resolving marker (global id marker).
Returns
ErrCode The error code indicating success or failure.

Definition at line 2334 of file iMOAB.cpp.

2335 { 2336  appData& data = context.appDatas[*pid]; 2337  ParallelComm* pco = context.pcomms[*pid]; 2338  EntityHandle cset = data.file_set; 2339  int dum_id = 0; 2340  ErrorCode rval; 2341  if( data.primary_elems.empty() ) 2342  { 2343  // skip actual resolve, assume vertices are distributed already , 2344  // no need to share them 2345  } 2346  else 2347  { 2348  // create an integer tag for resolving ; maybe it can be a long tag in the future 2349  // (more than 2 B vertices;) 2350  2351  Tag stag; 2352  rval = context.MBI->tag_get_handle( "__sharedmarker", 1, MB_TYPE_INTEGER, stag, MB_TAG_CREAT | MB_TAG_DENSE, 2353  &dum_id );MB_CHK_ERR( rval ); 2354  2355  if( *num_verts > (int)data.local_verts.size() ) 2356  { 2357  return moab::MB_FAILURE; 2358  } // we are not setting the size 2359  2360  rval = context.MBI->tag_set_data( stag, data.local_verts, (void*)marker );MB_CHK_ERR( rval ); // assumes integer tag 2361  2362  rval = pco->resolve_shared_ents( cset, -1, -1, &stag );MB_CHK_ERR( rval ); 2363  2364  rval = context.MBI->tag_delete( stag );MB_CHK_ERR( rval ); 2365  } 2366  // provide partition tag equal to rank 2367  Tag part_tag; 2368  dum_id = -1; 2369  rval = context.MBI->tag_get_handle( "PARALLEL_PARTITION", 1, MB_TYPE_INTEGER, part_tag, 2370  MB_TAG_CREAT | MB_TAG_SPARSE, &dum_id ); 2371  2372  if( part_tag == NULL || ( ( rval != MB_SUCCESS ) && ( rval != MB_ALREADY_ALLOCATED ) ) ) 2373  { 2374  std::cout << " can't get par part tag.\n"; 2375  return moab::MB_FAILURE; 2376  } 2377  2378  int rank = pco->rank(); 2379  rval = context.MBI->tag_set_data( part_tag, &cset, 1, &rank );MB_CHK_ERR( rval ); 2380  2381  return moab::MB_SUCCESS; 2382 }

References GlobalContext::appDatas, context, moab::Range::empty(), ErrorCode, appData::file_set, appData::local_verts, MB_ALREADY_ALLOCATED, MB_CHK_ERR, MB_SUCCESS, MB_TAG_CREAT, MB_TAG_DENSE, MB_TAG_SPARSE, MB_TYPE_INTEGER, GlobalContext::MBI, GlobalContext::pcomms, appData::primary_elems, moab::ParallelComm::rank(), moab::ParallelComm::resolve_shared_ents(), moab::Range::size(), moab::Interface::tag_delete(), moab::Interface::tag_get_handle(), and moab::Interface::tag_set_data().

◆ iMOAB_SetDoubleTagStorage()

ErrCode iMOAB_SetDoubleTagStorage ( iMOAB_AppID  pid,
const iMOAB_String  tag_storage_name,
int *  num_tag_storage_length,
int *  entity_type,
double *  tag_storage_data 
)

Store the specified values in a MOAB double Tag.

Note
Operations: Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]tag_storage_name(iMOAB_String) The tag names, separated, to store the data.
[in]num_tag_storage_length(int*) The size of total tag storage data (e.g., num_visible_vertices*components_per_entity or num_visible_elements*components_per_entity*num_tags).
[in]entity_type(int*) Type=0 for vertices, and Type=1 for primary elements.
[out]tag_storage_data(double*) The array data of type double to replace the internal tag memory; The data is assumed to be contiguous over the local set of visible entities (either vertices or elements). unrolled by tags
Returns
ErrCode The error code indicating success or failure.

Definition at line 1720 of file iMOAB.cpp.

1725 { 1726  ErrorCode rval; 1727  std::string tag_names( tag_storage_names ); 1728  // exactly the same code as for int tag :) maybe should check the type of tag too 1729  std::vector< std::string > tagNames; 1730  std::vector< Tag > tagHandles; 1731  std::string separator( ":" ); 1732  split_tag_names( tag_names, separator, tagNames ); 1733  1734  appData& data = context.appDatas[*pid]; 1735  Range* ents_to_set = NULL; 1736  1737  if( *ent_type == 0 ) // vertices 1738  { 1739  ents_to_set = &data.all_verts; 1740  } 1741  else if( *ent_type == 1 ) 1742  { 1743  ents_to_set = &data.primary_elems; 1744  } 1745  1746  int nents_to_be_set = (int)( *ents_to_set ).size(); 1747  int position = 0; 1748  1749  for( size_t i = 0; i < tagNames.size(); i++ ) 1750  { 1751  if( data.tagMap.find( tagNames[i] ) == data.tagMap.end() ) 1752  { 1753  return moab::MB_FAILURE; 1754  } // some tag not defined yet in the app 1755  1756  Tag tag = data.tagMap[tagNames[i]]; 1757  1758  int tagLength = 0; 1759  rval = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval ); 1760  1761  DataType dtype; 1762  rval = context.MBI->tag_get_data_type( tag, dtype );MB_CHK_ERR( rval ); 1763  1764  if( dtype != MB_TYPE_DOUBLE ) 1765  { 1766  return moab::MB_FAILURE; 1767  } 1768  1769  // set it on the subset of entities, based on type and length 1770  if( position + tagLength * nents_to_be_set > *num_tag_storage_length ) 1771  return moab::MB_FAILURE; // too many entity values to be set 1772  1773  rval = context.MBI->tag_set_data( tag, *ents_to_set, &tag_storage_data[position] );MB_CHK_ERR( rval ); 1774  // increment position to next tag 1775  position = position + tagLength * nents_to_be_set; 1776  } 1777  return moab::MB_SUCCESS; // no error 1778 }

References appData::all_verts, GlobalContext::appDatas, context, ErrorCode, MB_CHK_ERR, MB_SUCCESS, MB_TYPE_DOUBLE, GlobalContext::MBI, appData::primary_elems, moab::Range::size(), split_tag_names(), moab::Interface::tag_get_data_type(), moab::Interface::tag_get_length(), moab::Interface::tag_set_data(), and appData::tagMap.

◆ iMOAB_SetDoubleTagStorageWithGid()

ErrCode iMOAB_SetDoubleTagStorageWithGid ( iMOAB_AppID  pid,
const iMOAB_String  tag_storage_names,
int *  num_tag_storage_length,
int *  entity_type,
double *  tag_storage_data,
int *  globalIds 
)

Store the specified values in a MOAB double Tag, for given ids.

Note
Operations: Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]tag_storage_name(iMOAB_String) The tag names, separated, to store the data.
[in]num_tag_storage_length(int*) The size of total tag storage data (e.g., num_visible_vertices*components_per_entity or num_visible_elements*components_per_entity*num_tags).
[in]entity_type(int*) Type=0 for vertices, and Type=1 for primary elements.
[in]tag_storage_data(double*) The array data of type double to replace the internal tag memory; The data is assumed to be permuted over the local set of visible entities (either vertices or elements). unrolled by tags in parallel, might be different order, or different entities on the task
[in]globalIdsglobal ids of the cells to be set;
Returns
ErrCode The error code indicating success or failure.

Definition at line 1780 of file iMOAB.cpp.

1786 { 1787  ErrorCode rval; 1788  std::string tag_names( tag_storage_names ); 1789  // exactly the same code as for int tag :) maybe should check the type of tag too 1790  std::vector< std::string > tagNames; 1791  std::vector< Tag > tagHandles; 1792  std::string separator( ":" ); 1793  split_tag_names( tag_names, separator, tagNames ); 1794  1795  appData& data = context.appDatas[*pid]; 1796  Range* ents_to_set = NULL; 1797  1798  if( *ent_type == 0 ) // vertices 1799  { 1800  ents_to_set = &data.all_verts; 1801  } 1802  else if( *ent_type == 1 ) 1803  { 1804  ents_to_set = &data.primary_elems; 1805  } 1806  1807  int nents_to_be_set = (int)( *ents_to_set ).size(); 1808  1809  Tag gidTag = context.MBI->globalId_tag(); 1810  std::vector< int > gids; 1811  gids.resize( nents_to_be_set ); 1812  rval = context.MBI->tag_get_data( gidTag, *ents_to_set, &gids[0] );MB_CHK_ERR( rval ); 1813  1814  // so we will need to set the tags according to the global id passed; 1815  // so the order in tag_storage_data is the same as the order in globalIds, but the order 1816  // in local range is gids 1817  std::map< int, EntityHandle > eh_by_gid; 1818  int i = 0; 1819  for( Range::iterator it = ents_to_set->begin(); it != ents_to_set->end(); ++it, ++i ) 1820  { 1821  eh_by_gid[gids[i]] = *it; 1822  } 1823  1824  std::vector< int > tagLengths( tagNames.size() ); 1825  std::vector< Tag > tagList; 1826  size_t total_tag_len = 0; 1827  for( size_t i = 0; i < tagNames.size(); i++ ) 1828  { 1829  if( data.tagMap.find( tagNames[i] ) == data.tagMap.end() ) 1830  { 1831  MB_SET_ERR( moab::MB_FAILURE, "tag missing" ); 1832  } // some tag not defined yet in the app 1833  1834  Tag tag = data.tagMap[tagNames[i]]; 1835  tagList.push_back( tag ); 1836  1837  int tagLength = 0; 1838  rval = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval ); 1839  1840  total_tag_len += tagLength; 1841  tagLengths[i] = tagLength; 1842  DataType dtype; 1843  rval = context.MBI->tag_get_data_type( tag, dtype );MB_CHK_ERR( rval ); 1844  1845  if( dtype != MB_TYPE_DOUBLE ) 1846  { 1847  MB_SET_ERR( moab::MB_FAILURE, "tag not double type" ); 1848  } 1849  } 1850  bool serial = true; 1851 #ifdef MOAB_HAVE_MPI 1852  ParallelComm* pco = context.pcomms[*pid]; 1853  unsigned num_procs = pco->size(); 1854  if( num_procs > 1 ) serial = false; 1855 #endif 1856  1857  if( serial ) 1858  { 1859  // we do not assume anymore that the number of entities has to match 1860  // we will set only what matches, and skip entities that do not have corresponding global ids 1861  //assert( total_tag_len * nents_to_be_set - *num_tag_storage_length == 0 ); 1862  // tags are unrolled, we loop over global ids first, then careful about tags 1863  for( int i = 0; i < nents_to_be_set; i++ ) 1864  { 1865  int gid = globalIds[i]; 1866  std::map< int, EntityHandle >::iterator mapIt = eh_by_gid.find( gid ); 1867  if( mapIt == eh_by_gid.end() ) continue; 1868  EntityHandle eh = mapIt->second; 1869  // now loop over tags 1870  int indexInTagValues = 0; // 1871  for( size_t j = 0; j < tagList.size(); j++ ) 1872  { 1873  indexInTagValues += i * tagLengths[j]; 1874  rval = context.MBI->tag_set_data( tagList[j], &eh, 1, &tag_storage_data[indexInTagValues] );MB_CHK_ERR( rval ); 1875  // advance the pointer/index 1876  indexInTagValues += ( nents_to_be_set - i ) * tagLengths[j]; // at the end of tag data 1877  } 1878  } 1879  } 1880 #ifdef MOAB_HAVE_MPI 1881  else // it can be not serial only if pco->size() > 1, parallel 1882  { 1883  // in this case, we have to use 2 crystal routers, to send data to the processor that needs it 1884  // we will create first a tuple to rendevous points, then from there send to the processor that requested it 1885  // it is a 2-hop global gather scatter 1886  // TODO: allow for tags of different length; this is wrong 1887  int nbLocalVals = *num_tag_storage_length / ( (int)tagNames.size() ); // assumes all tags have the same length? 1888  // we do not expect the sizes to match 1889  //assert( nbLocalVals * tagNames.size() - *num_tag_storage_length == 0 ); 1890  TupleList TLsend; 1891  TLsend.initialize( 2, 0, 0, total_tag_len, nbLocalVals ); // to proc, marker(gid), total_tag_len doubles 1892  TLsend.enableWriteAccess(); 1893  // the processor id that processes global_id is global_id / num_ents_per_proc 1894  1895  int indexInRealLocal = 0; 1896  for( int i = 0; i < nbLocalVals; i++ ) 1897  { 1898  // to proc, marker, element local index, index in el 1899  int marker = globalIds[i]; 1900  int to_proc = marker % num_procs; 1901  int n = TLsend.get_n(); 1902  TLsend.vi_wr[2 * n] = to_proc; // send to processor 1903  TLsend.vi_wr[2 * n + 1] = marker; 1904  int indexInTagValues = 0; 1905  // tag data collect by number of tags 1906  for( size_t j = 0; j < tagList.size(); j++ ) 1907  { 1908  indexInTagValues += i * tagLengths[j]; 1909  for( int k = 0; k < tagLengths[j]; k++ ) 1910  { 1911  TLsend.vr_wr[indexInRealLocal++] = tag_storage_data[indexInTagValues + k]; 1912  } 1913  indexInTagValues += ( nbLocalVals - i ) * tagLengths[j]; 1914  } 1915  TLsend.inc_n(); 1916  } 1917  1918  //assert( nbLocalVals * total_tag_len - indexInRealLocal == 0 ); 1919  // send now requests, basically inform the rendez-vous point who needs a particular global id 1920  // send the data to the other processors: 1921  ( pco->proc_config().crystal_router() )->gs_transfer( 1, TLsend, 0 ); 1922  TupleList TLreq; 1923  TLreq.initialize( 2, 0, 0, 0, nents_to_be_set ); 1924  TLreq.enableWriteAccess(); 1925  for( int i = 0; i < nents_to_be_set; i++ ) 1926  { 1927  // to proc, marker 1928  int marker = gids[i]; 1929  int to_proc = marker % num_procs; 1930  int n = TLreq.get_n(); 1931  TLreq.vi_wr[2 * n] = to_proc; // send to processor 1932  TLreq.vi_wr[2 * n + 1] = marker; 1933  // tag data collect by number of tags 1934  TLreq.inc_n(); 1935  } 1936  ( pco->proc_config().crystal_router() )->gs_transfer( 1, TLreq, 0 ); 1937  1938  // we know now that process TLreq.vi_wr[2 * n] needs tags for gid TLreq.vi_wr[2 * n + 1] 1939  // we should first order by global id, and then build the new TL with send to proc, global id and 1940  // tags for it 1941  // sort by global ids the tuple lists 1942  moab::TupleList::buffer sort_buffer; 1943  sort_buffer.buffer_init( TLreq.get_n() ); 1944  TLreq.sort( 1, &sort_buffer ); 1945  sort_buffer.reset(); 1946  sort_buffer.buffer_init( TLsend.get_n() ); 1947  TLsend.sort( 1, &sort_buffer ); 1948  sort_buffer.reset(); 1949  // now send the tag values to the proc that requested it 1950  // in theory, for a full partition, TLreq and TLsend should have the same size, and 1951  // each dof should have exactly one target proc. Is that true or not in general ? 1952  // how do we plan to use this? Is it better to store the comm graph for future 1953  1954  // start copy from comm graph settle 1955  TupleList TLBack; 1956  TLBack.initialize( 3, 0, 0, total_tag_len, 0 ); // to proc, marker, tag from proc , tag values 1957  TLBack.enableWriteAccess(); 1958  1959  int n1 = TLreq.get_n(); 1960  int n2 = TLsend.get_n(); 1961  1962  int indexInTLreq = 0; 1963  int indexInTLsend = 0; // advance both, according to the marker 1964  if( n1 > 0 && n2 > 0 ) 1965  { 1966  1967  while( indexInTLreq < n1 && indexInTLsend < n2 ) // if any is over, we are done 1968  { 1969  int currentValue1 = TLreq.vi_rd[2 * indexInTLreq + 1]; 1970  int currentValue2 = TLsend.vi_rd[2 * indexInTLsend + 1]; 1971  if( currentValue1 < currentValue2 ) 1972  { 1973  // we have a big problem; basically, we are saying that 1974  // dof currentValue is on one model and not on the other 1975  // std::cout << " currentValue1:" << currentValue1 << " missing in comp2" << "\n"; 1976  indexInTLreq++; 1977  continue; 1978  } 1979  if( currentValue1 > currentValue2 ) 1980  { 1981  // std::cout << " currentValue2:" << currentValue2 << " missing in comp1" << "\n"; 1982  indexInTLsend++; 1983  continue; 1984  } 1985  int size1 = 1; 1986  int size2 = 1; 1987  while( indexInTLreq + size1 < n1 && currentValue1 == TLreq.vi_rd[2 * ( indexInTLreq + size1 ) + 1] ) 1988  size1++; 1989  while( indexInTLsend + size2 < n2 && currentValue2 == TLsend.vi_rd[2 * ( indexInTLsend + size2 ) + 1] ) 1990  size2++; 1991  // must be found in both lists, find the start and end indices 1992  for( int i1 = 0; i1 < size1; i1++ ) 1993  { 1994  for( int i2 = 0; i2 < size2; i2++ ) 1995  { 1996  // send the info back to components 1997  int n = TLBack.get_n(); 1998  TLBack.reserve(); 1999  TLBack.vi_wr[3 * n] = TLreq.vi_rd[2 * ( indexInTLreq + i1 )]; // send back to the proc marker 2000  // came from, info from comp2 2001  TLBack.vi_wr[3 * n + 1] = currentValue1; // initial value (resend, just for verif ?) 2002  TLBack.vi_wr[3 * n + 2] = TLsend.vi_rd[2 * ( indexInTLsend + i2 )]; // from proc on comp2 2003  // also fill tag values 2004  for( size_t k = 0; k < total_tag_len; k++ ) 2005  { 2006  TLBack.vr_rd[total_tag_len * n + k] = 2007  TLsend.vr_rd[total_tag_len * indexInTLsend + k]; // deep copy of tag values 2008  } 2009  } 2010  } 2011  indexInTLreq += size1; 2012  indexInTLsend += size2; 2013  } 2014  } 2015  ( pco->proc_config().crystal_router() )->gs_transfer( 1, TLBack, 0 ); 2016  // end copy from comm graph 2017  // after we are done sending, we need to set those tag values, in a reverse process compared to send 2018  n1 = TLBack.get_n(); 2019  double* ptrVal = &TLBack.vr_rd[0]; // 2020  for( int i = 0; i < n1; i++ ) 2021  { 2022  int gid = TLBack.vi_rd[3 * i + 1]; // marker 2023  std::map< int, EntityHandle >::iterator mapIt = eh_by_gid.find( gid ); 2024  if( mapIt == eh_by_gid.end() ) continue; 2025  EntityHandle eh = mapIt->second; 2026  // now loop over tags 2027  2028  for( size_t j = 0; j < tagList.size(); j++ ) 2029  { 2030  rval = context.MBI->tag_set_data( tagList[j], &eh, 1, (void*)ptrVal );MB_CHK_ERR( rval ); 2031  // advance the pointer/index 2032  ptrVal += tagLengths[j]; // at the end of tag data per call 2033  } 2034  } 2035  } 2036 #endif 2037  return MB_SUCCESS; 2038 }

References appData::all_verts, GlobalContext::appDatas, moab::Range::begin(), context, moab::ProcConfig::crystal_router(), moab::TupleList::enableWriteAccess(), moab::Range::end(), ErrorCode, moab::TupleList::get_n(), moab::Interface::globalId_tag(), moab::TupleList::inc_n(), moab::TupleList::initialize(), MB_CHK_ERR, MB_SET_ERR, MB_SUCCESS, MB_TYPE_DOUBLE, GlobalContext::MBI, GlobalContext::pcomms, appData::primary_elems, moab::ParallelComm::proc_config(), moab::TupleList::reserve(), moab::TupleList::buffer::reset(), moab::Range::size(), moab::ParallelComm::size(), moab::TupleList::sort(), split_tag_names(), moab::Interface::tag_get_data(), moab::Interface::tag_get_data_type(), moab::Interface::tag_get_length(), moab::Interface::tag_set_data(), appData::tagMap, moab::TupleList::vi_rd, moab::TupleList::vi_wr, moab::TupleList::vr_rd, and moab::TupleList::vr_wr.

◆ iMOAB_SetGlobalInfo()

ErrCode iMOAB_SetGlobalInfo ( iMOAB_AppID  pid,
int *  num_global_verts,
int *  num_global_elems 
)

Set global information for number of vertices and number of elements; It is usually available from the h5m file, or it can be computed with a MPI_Reduce.

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]num_global_verts(int*) The number of total vertices in the mesh.
[in]num_global_elems(MPI_Comm) The number of total elements in the mesh.
Returns
ErrCode The error code indicating success or failure.

Definition at line 2308 of file iMOAB.cpp.

2309 { 2310  appData& data = context.appDatas[*pid]; 2311  data.num_global_vertices = *num_global_verts; 2312  data.num_global_elements = *num_global_elems; 2313  return moab::MB_SUCCESS; 2314 }

References GlobalContext::appDatas, context, MB_SUCCESS, appData::num_global_elements, and appData::num_global_vertices.

Referenced by iMOAB_UpdateMeshInfo().

◆ iMOAB_SetIntTagStorage()

ErrCode iMOAB_SetIntTagStorage ( iMOAB_AppID  pid,
const iMOAB_String  tag_storage_name,
int *  num_tag_storage_length,
int *  entity_type,
int *  tag_storage_data 
)

Store the specified values in a MOAB integer Tag.

Values are set on vertices or elements, depending on entity_type

Note
we could change the api to accept as input the tag index, as returned by iMOAB_DefineTagStorage.

Operations: Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]tag_storage_name(iMOAB_String) The tag name to store/retreive the data in MOAB.
[in]num_tag_storage_length(int*) The size of tag storage data (e.g., num_visible_vertices*components_per_entity or num_visible_elements*components_per_entity).
[in]entity_type(int*) Type=0 for vertices, and Type=1 for primary elements.
[out]tag_storage_data(int*) The array data of type int to replace the internal tag memory; The data is assumed to be contiguous over the local set of visible entities (either vertices or elements).
Returns
ErrCode The error code indicating success or failure.

Definition at line 1606 of file iMOAB.cpp.

1611 { 1612  std::string tag_name( tag_storage_name ); 1613  1614  // Get the application data 1615  appData& data = context.appDatas[*pid]; 1616  1617  if( data.tagMap.find( tag_name ) == data.tagMap.end() ) 1618  { 1619  return moab::MB_FAILURE; 1620  } // tag not defined 1621  1622  Tag tag = data.tagMap[tag_name]; 1623  1624  int tagLength = 0; 1625  ErrorCode rval = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval ); 1626  1627  DataType dtype; 1628  rval = context.MBI->tag_get_data_type( tag, dtype ); 1629  1630  if( MB_SUCCESS != rval || dtype != MB_TYPE_INTEGER ) 1631  { 1632  return moab::MB_FAILURE; 1633  } 1634  1635  // set it on a subset of entities, based on type and length 1636  Range* ents_to_set; 1637  1638  if( *ent_type == 0 ) // vertices 1639  { 1640  ents_to_set = &data.all_verts; 1641  } 1642  else // if (*ent_type == 1) // *ent_type can be 0 (vertices) or 1 (elements) 1643  { 1644  ents_to_set = &data.primary_elems; 1645  } 1646  1647  int nents_to_be_set = *num_tag_storage_length / tagLength; 1648  1649  if( nents_to_be_set > (int)ents_to_set->size() ) 1650  { 1651  return moab::MB_FAILURE; 1652  } // to many entities to be set or too few 1653  1654  // restrict the range; everything is contiguous; or not? 1655  1656  // Range contig_range ( * ( ents_to_set->begin() ), * ( ents_to_set->begin() + nents_to_be_set - 1657  // 1 ) ); 1658  rval = context.MBI->tag_set_data( tag, *ents_to_set, tag_storage_data );MB_CHK_ERR( rval ); 1659  1660  return moab::MB_SUCCESS; // no error 1661 }

References appData::all_verts, GlobalContext::appDatas, context, ErrorCode, MB_CHK_ERR, MB_SUCCESS, MB_TYPE_INTEGER, GlobalContext::MBI, appData::primary_elems, moab::Range::size(), moab::Interface::tag_get_data_type(), moab::Interface::tag_get_length(), moab::Interface::tag_set_data(), and appData::tagMap.

◆ iMOAB_SynchronizeTags()

ErrCode iMOAB_SynchronizeTags ( iMOAB_AppID  pid,
int *  num_tag,
int *  tag_indices,
int *  ent_type 
)

Exchange tag values for the given tags across process boundaries.

Note
Operations: Collective
Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]num_tags(int*) Number of tags to exchange.
[in]tag_indices(int*) Array with tag indices of interest (size = *num_tags).
[in]ent_type(int*) The type of entity for tag exchange.
Returns
ErrCode The error code indicating success or failure.

Definition at line 2099 of file iMOAB.cpp.

2100 { 2101 #ifdef MOAB_HAVE_MPI 2102  appData& data = context.appDatas[*pid]; 2103  Range ent_exchange; 2104  std::vector< Tag > tags; 2105  2106  for( int i = 0; i < *num_tag; i++ ) 2107  { 2108  if( tag_indices[i] < 0 || tag_indices[i] >= (int)data.tagList.size() ) 2109  { 2110  return moab::MB_FAILURE; 2111  } // error in tag index 2112  2113  tags.push_back( data.tagList[tag_indices[i]] ); 2114  } 2115  2116  if( *ent_type == 0 ) 2117  { 2118  ent_exchange = data.all_verts; 2119  } 2120  else if( *ent_type == 1 ) 2121  { 2122  ent_exchange = data.primary_elems; 2123  } 2124  else 2125  { 2126  return moab::MB_FAILURE; 2127  } // unexpected type 2128  2129  ParallelComm* pco = context.pcomms[*pid]; 2130  2131  ErrorCode rval = pco->exchange_tags( tags, tags, ent_exchange );MB_CHK_ERR( rval ); 2132  2133 #else 2134  /* do nothing if serial */ 2135  // just silence the warning 2136  // do not call sync tags in serial! 2137  int k = *pid + *num_tag + *tag_indices + *ent_type; 2138  k++; 2139 #endif 2140  2141  return moab::MB_SUCCESS; 2142 }

References appData::all_verts, GlobalContext::appDatas, context, ErrorCode, moab::ParallelComm::exchange_tags(), MB_CHK_ERR, MB_SUCCESS, GlobalContext::pcomms, appData::primary_elems, and appData::tagList.

◆ iMOAB_UpdateMeshInfo()

ErrCode iMOAB_UpdateMeshInfo ( iMOAB_AppID  pid)

Update local mesh data structure, from file information.

Note
The method should be called after mesh modifications, for example reading a file or creating mesh in memory

Operations: Collective

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
Returns
ErrCode The error code indicating success or failure.

Definition at line 782 of file iMOAB.cpp.

783 { 784  // this will include ghost elements info 785  appData& data = context.appDatas[*pid]; 786  EntityHandle fileSet = data.file_set; 787  788  // first clear all data ranges; this can be called after ghosting 789  data.all_verts.clear(); 790  data.primary_elems.clear(); 791  data.local_verts.clear(); 792  data.owned_verts.clear(); 793  data.ghost_vertices.clear(); 794  data.owned_elems.clear(); 795  data.ghost_elems.clear(); 796  data.mat_sets.clear(); 797  data.neu_sets.clear(); 798  data.diri_sets.clear(); 799  800  // Let us get all the vertex entities 801  ErrorCode rval = context.MBI->get_entities_by_type( fileSet, MBVERTEX, data.all_verts, true );MB_CHK_ERR( rval ); // recursive 802  803  // Let us check first entities of dimension = 3 804  data.dimension = 3; 805  rval = context.MBI->get_entities_by_dimension( fileSet, data.dimension, data.primary_elems, true );MB_CHK_ERR( rval ); // recursive 806  807  if( data.primary_elems.empty() ) 808  { 809  // Now 3-D elements. Let us check entities of dimension = 2 810  data.dimension = 2; 811  rval = context.MBI->get_entities_by_dimension( fileSet, data.dimension, data.primary_elems, true );MB_CHK_ERR( rval ); // recursive 812  813  if( data.primary_elems.empty() ) 814  { 815  // Now 3-D/2-D elements. Let us check entities of dimension = 1 816  data.dimension = 1; 817  rval = context.MBI->get_entities_by_dimension( fileSet, data.dimension, data.primary_elems, true );MB_CHK_ERR( rval ); // recursive 818  819  if( data.primary_elems.empty() ) 820  { 821  // no elements of dimension 1 or 2 or 3; it could happen for point clouds 822  data.dimension = 0; 823  } 824  } 825  } 826  827  // check if the current mesh is just a point cloud 828  data.point_cloud = ( ( data.primary_elems.size() == 0 && data.all_verts.size() > 0 ) || data.dimension == 0 ); 829  830 #ifdef MOAB_HAVE_MPI 831  832  if( context.MPI_initialized ) 833  { 834  ParallelComm* pco = context.pcomms[*pid]; 835  836  // filter ghost vertices, from local 837  rval = pco->filter_pstatus( data.all_verts, PSTATUS_GHOST, PSTATUS_NOT, -1, &data.local_verts );MB_CHK_ERR( rval ); 838  839  // Store handles for all ghosted entities 840  data.ghost_vertices = subtract( data.all_verts, data.local_verts ); 841  842  // filter ghost elements, from local 843  rval = pco->filter_pstatus( data.primary_elems, PSTATUS_GHOST, PSTATUS_NOT, -1, &data.owned_elems );MB_CHK_ERR( rval ); 844  845  data.ghost_elems = subtract( data.primary_elems, data.owned_elems ); 846  // now update global number of primary cells and global number of vertices 847  // determine first number of owned vertices 848  // Get local owned vertices 849  rval = pco->filter_pstatus( data.all_verts, PSTATUS_NOT_OWNED, PSTATUS_NOT, -1, &data.owned_verts );MB_CHK_ERR( rval ); 850  int local[2], global[2]; 851  local[0] = data.owned_verts.size(); 852  local[1] = data.owned_elems.size(); 853  MPI_Allreduce( local, global, 2, MPI_INT, MPI_SUM, pco->comm() ); 854  rval = iMOAB_SetGlobalInfo( pid, &( global[0] ), &( global[1] ) );MB_CHK_ERR( rval ); 855  } 856  else 857  { 858  data.local_verts = data.all_verts; 859  data.owned_elems = data.primary_elems; 860  } 861  862 #else 863  864  data.local_verts = data.all_verts; 865  data.owned_elems = data.primary_elems; 866  867 #endif 868  869  // Get the references for some standard internal tags such as material blocks, BCs, etc 870  rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.material_tag ), 0, 1, 871  data.mat_sets, Interface::UNION );MB_CHK_ERR( rval ); 872  873  rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.neumann_tag ), 0, 1, 874  data.neu_sets, Interface::UNION );MB_CHK_ERR( rval ); 875  876  rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.dirichlet_tag ), 0, 1, 877  data.diri_sets, Interface::UNION );MB_CHK_ERR( rval ); 878  879  return moab::MB_SUCCESS; 880 }

References appData::all_verts, GlobalContext::appDatas, moab::Range::clear(), moab::ParallelComm::comm(), context, appData::dimension, appData::diri_sets, GlobalContext::dirichlet_tag, moab::Range::empty(), ErrorCode, appData::file_set, moab::ParallelComm::filter_pstatus(), moab::Interface::get_entities_by_dimension(), moab::Interface::get_entities_by_type(), moab::Interface::get_entities_by_type_and_tag(), appData::ghost_elems, appData::ghost_vertices, iMOAB_SetGlobalInfo(), appData::local_verts, appData::mat_sets, GlobalContext::material_tag, MB_CHK_ERR, MB_SUCCESS, MBENTITYSET, GlobalContext::MBI, MBVERTEX, GlobalContext::MPI_initialized, appData::neu_sets, GlobalContext::neumann_tag, appData::owned_elems, appData::owned_verts, GlobalContext::pcomms, appData::point_cloud, appData::primary_elems, PSTATUS_GHOST, PSTATUS_NOT, PSTATUS_NOT_OWNED, moab::Range::size(), moab::subtract(), and moab::Interface::UNION.

Referenced by iMOAB_DetermineGhostEntities(), iMOAB_LoadMesh(), iMOAB_MergeVertices(), and iMOAB_ReceiveMesh().

◆ iMOAB_WriteLocalMesh()

ErrCode iMOAB_WriteLocalMesh ( iMOAB_AppID  pid,
iMOAB_String  prefix 
)

Write a local MOAB mesh copy.

Note
The interface will write one file per task. Only the local mesh set and solution data associated to the application will be written to the file. Very convenient method used for debugging.

Operations: Not Collective (independent operation).

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]prefix(iMOAB_String) The MOAB file prefix that will be used with task ID in parallel.
Returns
ErrCode The error code indicating success or failure.

Definition at line 764 of file iMOAB.cpp.

765 { 766  IMOAB_CHECKPOINTER( prefix, 2 ); 767  IMOAB_ASSERT( strlen( prefix ), "Invalid prefix string length." ); 768  769  std::ostringstream file_name; 770  int rank = 0, size = 1; 771 #ifdef MOAB_HAVE_MPI 772  rank = context.pcomms[*pid]->rank(); 773  size = context.pcomms[*pid]->size(); 774 #endif 775  file_name << prefix << "_" << size << "_" << rank << ".h5m"; 776  // Now let us actually write the file to disk with appropriate options 777  ErrorCode rval = context.MBI->write_file( file_name.str().c_str(), 0, 0, &context.appDatas[*pid].file_set, 1 );MB_CHK_ERR( rval ); 778  779  return moab::MB_SUCCESS; 780 }

References GlobalContext::appDatas, context, ErrorCode, IMOAB_ASSERT, IMOAB_CHECKPOINTER, MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, GlobalContext::pcomms, size, and moab::Interface::write_file().

◆ iMOAB_WriteMesh()

ErrCode iMOAB_WriteMesh ( iMOAB_AppID  pid,
const iMOAB_String  filename,
const iMOAB_String  write_options 
)

Write a MOAB mesh along with the solution tags to a file.

Note
The interface will write one single file (H5M) and for serial files (VTK/Exodus), it will write one file per task. Write options include parallel write options, if needed. Only the mesh set and solution data associated to the application will be written to the file.

Operations: Collective for parallel write, non collective for serial write.

Parameters
[in]pid(iMOAB_AppID) The unique pointer to the application ID.
[in]filename(iMOAB_String) The MOAB mesh file (H5M) to write all the entities contained in the internal application mesh set.
[in]write_options(iMOAB_String) Additional options for writing the MOAB mesh in parallel.
Returns
ErrCode The error code indicating success or failure.

Definition at line 703 of file iMOAB.cpp.

704 { 705  IMOAB_CHECKPOINTER( filename, 2 ); 706  IMOAB_ASSERT( strlen( filename ), "Invalid filename length." ); 707  708  appData& data = context.appDatas[*pid]; 709  EntityHandle fileSet = data.file_set; 710  711  std::ostringstream newopts; 712 #ifdef MOAB_HAVE_MPI 713  std::string write_opts( ( write_options ? write_options : "" ) ); 714  std::string pcid( "PARALLEL_COMM=" ); 715  std::size_t found = write_opts.find( pcid ); 716  717  if( found != std::string::npos ) 718  { 719  std::cerr << " cannot specify PARALLEL_COMM option, it is implicit \n"; 720  return moab::MB_FAILURE; 721  } 722  723  // if write in parallel, add pc option, to be sure about which ParallelComm instance is used 724  std::string pw( "PARALLEL=WRITE_PART" ); 725  found = write_opts.find( pw ); 726  727  if( found != std::string::npos ) 728  { 729  newopts << "PARALLEL_COMM=" << *pid << ";"; 730  } 731  732 #endif 733  734  if( write_options ) newopts << write_options; 735  736  std::vector< Tag > copyTagList = data.tagList; 737  // append Global ID and Parallel Partition 738  std::string gid_name_tag( "GLOBAL_ID" ); 739  740  // export global id tag, we need it always 741  if( data.tagMap.find( gid_name_tag ) == data.tagMap.end() ) 742  { 743  Tag gid = context.MBI->globalId_tag(); 744  copyTagList.push_back( gid ); 745  } 746  // also Parallel_Partition PARALLEL_PARTITION 747  std::string pp_name_tag( "PARALLEL_PARTITION" ); 748  749  // write parallel part tag too, if it exists 750  if( data.tagMap.find( pp_name_tag ) == data.tagMap.end() ) 751  { 752  Tag ptag; 753  context.MBI->tag_get_handle( pp_name_tag.c_str(), ptag ); 754  if( ptag ) copyTagList.push_back( ptag ); 755  } 756  757  // Now let us actually write the file to disk with appropriate options 758  ErrorCode rval = context.MBI->write_file( filename, 0, newopts.str().c_str(), &fileSet, 1, &copyTagList[0], 759  (int)copyTagList.size() );MB_CHK_ERR( rval ); 760  761  return moab::MB_SUCCESS; 762 }

References GlobalContext::appDatas, context, ErrorCode, appData::file_set, moab::Interface::globalId_tag(), IMOAB_ASSERT, IMOAB_CHECKPOINTER, MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, moab::Interface::tag_get_handle(), appData::tagList, appData::tagMap, and moab::Interface::write_file().