MOAB: Mesh Oriented datABase  (version 5.4.1)
iMOAB.cpp File Reference
#include "moab/MOABConfig.h"
#include "moab/Core.hpp"
#include "DebugOutput.hpp"
#include "moab/iMOAB.h"
#include "moab/CartVect.hpp"
#include "MBTagConventions.hpp"
#include "moab/MeshTopoUtil.hpp"
#include "moab/ReadUtilIface.hpp"
#include "moab/MergeMesh.hpp"
#include <cassert>
#include <sstream>
#include <iostream>
+ Include dependency graph for iMOAB.cpp:

Go to the source code of this file.

Classes

struct  appData
struct  GlobalContext

Functions

ErrCode iMOAB_Initialize (int argc, iMOAB_String *argv)
 Initialize the iMOAB interface implementation.
ErrCode iMOAB_InitializeFortran ()
 Initialize the iMOAB interface implementation from Fortran driver.
ErrCode iMOAB_Finalize ()
 Finalize the iMOAB interface implementation.
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.
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.
ErrCode iMOAB_DeregisterApplication (iMOAB_AppID pid)
 De-Register application: delete mesh (set) associated with the application ID.
ErrCode iMOAB_DeregisterApplicationFortran (iMOAB_AppID pid)
 De-Register the Fortran based application: delete mesh (set) associated with the application ID.
static void split_tag_names (std::string input_names, std::string &separator, std::vector< std::string > &list_tag_names)
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.
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.
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.
ErrCode iMOAB_WriteLocalMesh (iMOAB_AppID pid, iMOAB_String prefix)
 Write a local MOAB mesh copy.
ErrCode iMOAB_UpdateMeshInfo (iMOAB_AppID pid)
 Update local mesh data structure, from file information.
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.
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.
ErrCode iMOAB_GetVertexOwnership (iMOAB_AppID pid, int *vertices_length, int *visible_global_rank_ID)
 Get vertex ownership information.
ErrCode iMOAB_GetVisibleVerticesCoordinates (iMOAB_AppID pid, int *coords_length, double *coordinates)
 Get vertex coordinates for all local (owned and ghosted) vertices.
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.
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).
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.
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.
ErrCode iMOAB_GetElementConnectivity (iMOAB_AppID pid, iMOAB_LocalID *elem_index, int *connectivity_length, int *element_connectivity)
 Get the connectivity for one element only.
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.
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.
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.
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.
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.
ErrCode iMOAB_SetIntTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_name, int *num_tag_storage_length, int *ent_type, int *tag_storage_data)
 Store the specified values in a MOAB integer Tag.
ErrCode iMOAB_GetIntTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_name, int *num_tag_storage_length, int *ent_type, int *tag_storage_data)
 Get the specified values in a MOAB integer Tag.
ErrCode iMOAB_SetDoubleTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_names, int *num_tag_storage_length, int *ent_type, double *tag_storage_data)
 Store the specified values in a MOAB double Tag.
ErrCode iMOAB_SetDoubleTagStorageWithGid (iMOAB_AppID pid, const iMOAB_String tag_storage_names, int *num_tag_storage_length, int *ent_type, double *tag_storage_data, int *globalIds)
 Store the specified values in a MOAB double Tag, for given ids.
ErrCode iMOAB_GetDoubleTagStorage (iMOAB_AppID pid, const iMOAB_String tag_storage_names, int *num_tag_storage_length, int *ent_type, double *tag_storage_data)
 Retrieve the specified values in a MOAB double Tag.
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.
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.
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.
ErrCode iMOAB_CreateVertices (iMOAB_AppID pid, int *coords_len, int *dim, double *coordinates)
 Create vertices for an app; it assumes no other vertices.
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.
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.
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.

Variables

static struct GlobalContext context

Detailed Description

Definition in file iMOAB.cpp.


Function Documentation

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 2211 of file iMOAB.cpp.

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().

Referenced by main().

{
    // Create elements
    appData& data = context.appDatas[*pid];

    ReadUtilIface* read_iface;
    ErrorCode rval = context.MBI->query_interface( read_iface );MB_CHK_ERR( rval );

    EntityType mbtype = (EntityType)( *type );
    EntityHandle actual_start_handle;
    EntityHandle* array = NULL;
    rval = read_iface->get_element_connect( *num_elem, *num_nodes_per_element, mbtype, 1, actual_start_handle, array );MB_CHK_ERR( rval );

    // fill up with actual connectivity from input; assume the vertices are in order, and start
    // vertex is the first in the current data vertex range
    EntityHandle firstVertex = data.local_verts[0];

    for( int j = 0; j < *num_elem * ( *num_nodes_per_element ); j++ )
    {
        array[j] = connectivity[j] + firstVertex - 1;
    }  // assumes connectivity uses 1 based array (from fortran, mostly)

    Range new_elems( actual_start_handle, actual_start_handle + *num_elem - 1 );

    rval = context.MBI->add_entities( data.file_set, new_elems );MB_CHK_ERR( rval );

    data.primary_elems.merge( new_elems );

    // add to adjacency
    rval = read_iface->update_adjacencies( actual_start_handle, *num_elem, *num_nodes_per_element, array );MB_CHK_ERR( rval );
    // organize all new elements in block, with the given block ID; if the block set is not
    // existing, create  a new mesh set;
    Range sets;
    int set_no            = *block_ID;
    const void* setno_ptr = &set_no;
    rval = context.MBI->get_entities_by_type_and_tag( data.file_set, MBENTITYSET, &context.material_tag, &setno_ptr, 1,
                                                      sets );
    EntityHandle block_set;

    if( MB_FAILURE == rval || sets.empty() )
    {
        // create a new set, with this block ID
        rval = context.MBI->create_meshset( MESHSET_SET, block_set );MB_CHK_ERR( rval );

        rval = context.MBI->tag_set_data( context.material_tag, &block_set, 1, &set_no );MB_CHK_ERR( rval );

        // add the material set to file set
        rval = context.MBI->add_entities( data.file_set, &block_set, 1 );MB_CHK_ERR( rval );
    }
    else
    {
        block_set = sets[0];
    }  // first set is the one we want

    /// add the new ents to the clock set
    rval = context.MBI->add_entities( block_set, new_elems );MB_CHK_ERR( rval );

    return moab::MB_SUCCESS;
}
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 2190 of file iMOAB.cpp.

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().

Referenced by main().

{
    ErrorCode rval;
    appData& data = context.appDatas[*pid];

    if( !data.local_verts.empty() )  // we should have no vertices in the app
    {
        return moab::MB_FAILURE;
    }

    int nverts = *coords_len / *dim;

    rval = context.MBI->create_vertices( coordinates, nverts, data.local_verts );MB_CHK_ERR( rval );

    rval = context.MBI->add_entities( data.file_set, data.local_verts );MB_CHK_ERR( rval );

    // also add the vertices to the all_verts range
    data.all_verts.merge( data.local_verts );
    return moab::MB_SUCCESS;
}
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 1448 of file iMOAB.cpp.

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, moab::ParallelComm::rank(), split_tag_names(), moab::Interface::tag_get_handle(), appData::tagList, appData::tagMap, and TagType.

Referenced by commgraphtest(), main(), and migrate().

{
    // see if the tag is already existing, and if yes, check the type, length
    if( *tag_type < 0 || *tag_type > 5 )
    {
        return moab::MB_FAILURE;
    }  // we have 6 types of tags supported so far

    DataType tagDataType;
    TagType tagType;
    void* defaultVal        = NULL;
    int* defInt             = new int[*components_per_entity];
    double* defDouble       = new double[*components_per_entity];
    EntityHandle* defHandle = new EntityHandle[*components_per_entity];

    for( int i = 0; i < *components_per_entity; i++ )
    {
        defInt[i]    = 0;
        defDouble[i] = -1e+10;
        defHandle[i] = (EntityHandle)0;
    }

    switch( *tag_type )
    {
        case 0:
            tagDataType = MB_TYPE_INTEGER;
            tagType     = MB_TAG_DENSE;
            defaultVal  = defInt;
            break;

        case 1:
            tagDataType = MB_TYPE_DOUBLE;
            tagType     = MB_TAG_DENSE;
            defaultVal  = defDouble;
            break;

        case 2:
            tagDataType = MB_TYPE_HANDLE;
            tagType     = MB_TAG_DENSE;
            defaultVal  = defHandle;
            break;

        case 3:
            tagDataType = MB_TYPE_INTEGER;
            tagType     = MB_TAG_SPARSE;
            defaultVal  = defInt;
            break;

        case 4:
            tagDataType = MB_TYPE_DOUBLE;
            tagType     = MB_TAG_SPARSE;
            defaultVal  = defDouble;
            break;

        case 5:
            tagDataType = MB_TYPE_HANDLE;
            tagType     = MB_TAG_SPARSE;
            defaultVal  = defHandle;
            break;

        default: {
            delete[] defInt;
            delete[] defDouble;
            delete[] defHandle;
            return moab::MB_FAILURE;
        }  // error
    }

    Tag tagHandle;
    // split storage names if separated list

    std::string tag_name( tag_storage_name );

    //  first separate the names of the tags
    // we assume that there are separators ":" between the tag names
    std::vector< std::string > tagNames;
    std::string separator( ":" );
    split_tag_names( tag_name, separator, tagNames );

    ErrorCode rval = moab::MB_SUCCESS;  // assume success already :)
    appData& data  = context.appDatas[*pid];
    int already_defined_tags = 0;

    for( size_t i = 0; i < tagNames.size(); i++ )
    {
        rval = context.MBI->tag_get_handle( tagNames[i].c_str(), *components_per_entity, tagDataType, tagHandle,
                                            tagType | MB_TAG_EXCL | MB_TAG_CREAT, defaultVal );

        if( MB_ALREADY_ALLOCATED == rval )
        {
            std::map< std::string, Tag >& mTags        = data.tagMap;
            std::map< std::string, Tag >::iterator mit = mTags.find( tagNames[i].c_str() );

            if( mit == mTags.end() )
            {
                // add it to the map
                mTags[tagNames[i]] = tagHandle;
                // push it to the list of tags, too
                *tag_index = (int)data.tagList.size();
                data.tagList.push_back( tagHandle );
            }
            rval = MB_SUCCESS;
            already_defined_tags++;
        }
        else if( MB_SUCCESS == rval )
        {
            data.tagMap[tagNames[i]] = tagHandle;
            *tag_index               = (int)data.tagList.size();
            data.tagList.push_back( tagHandle );
        }
        else
        {
            rval = moab::MB_FAILURE;  // some tags were not created
        }
    }
    // we don't need default values anymore, avoid leaks
    int rankHere  = 0;
#ifdef MOAB_HAVE_MPI
    ParallelComm* pco = context.pcomms[*pid];
    rankHere          = pco->rank();
#endif
    if( !rankHere && already_defined_tags)
        std::cout << " application with ID: " << *pid << " global id: " << data.global_id << " name: " << data.name
                  << " has "  << already_defined_tags << " already defined tags out of " << tagNames.size() << " tags \n";
    delete[] defInt;
    delete[] defDouble;
    delete[] defHandle;
    return rval;
}

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 327 of file iMOAB.cpp.

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, moab::ParallelComm::rank(), moab::Range::subset_by_type(), moab::subtract(), moab::Interface::UNION, and GlobalContext::unused_pid.

Referenced by commgraphtest(), iMOAB_DeregisterApplicationFortran(), main(), migrate(), and migrate_smart().

{
    // the file set , parallel comm are all in vectors indexed by *pid
    // assume we did not delete anything yet
    // *pid will not be reused if we register another application
    appData& data = context.appDatas[*pid];
    int rankHere  = 0;
#ifdef MOAB_HAVE_MPI
    ParallelComm* pco = context.pcomms[*pid];
    rankHere          = pco->rank();
#endif
    if( !rankHere )
        std::cout << " application with ID: " << *pid << " global id: " << data.global_id << " name: " << data.name
                  << " is de-registered now \n";

    EntityHandle fileSet = data.file_set;
    // get all entities part of the file set
    Range fileents;
    ErrorCode rval = context.MBI->get_entities_by_handle( fileSet, fileents, /*recursive */ true );MB_CHK_ERR( rval );

    fileents.insert( fileSet );

    rval = context.MBI->get_entities_by_type( fileSet, MBENTITYSET, fileents );MB_CHK_ERR( rval );  // append all mesh sets

#ifdef MOAB_HAVE_TEMPESTREMAP
    if( data.tempestData.remapper ) delete data.tempestData.remapper;
    if( data.tempestData.weightMaps.size() ) data.tempestData.weightMaps.clear();
#endif

#ifdef MOAB_HAVE_MPI

    // we could get the pco also with
    // ParallelComm * pcomm = ParallelComm::get_pcomm(context.MBI, *pid);

    std::map< int, ParCommGraph* >& pargs = data.pgraph;

    // free the parallel comm graphs associated with this app
    for( std::map< int, ParCommGraph* >::iterator mt = pargs.begin(); mt != pargs.end(); mt++ )
    {
        ParCommGraph* pgr = mt->second;
        delete pgr;
        pgr = NULL;
    }
    if( pco ) delete pco;
#endif

    // delete first all except vertices
    Range vertices = fileents.subset_by_type( MBVERTEX );
    Range noverts  = subtract( fileents, vertices );

    rval = context.MBI->delete_entities( noverts );MB_CHK_ERR( rval );
    // now retrieve connected elements that still exist (maybe in other sets, pids?)
    Range adj_ents_left;
    rval = context.MBI->get_adjacencies( vertices, 1, false, adj_ents_left, Interface::UNION );MB_CHK_ERR( rval );
    rval = context.MBI->get_adjacencies( vertices, 2, false, adj_ents_left, Interface::UNION );MB_CHK_ERR( rval );
    rval = context.MBI->get_adjacencies( vertices, 3, false, adj_ents_left, Interface::UNION );MB_CHK_ERR( rval );

    if( !adj_ents_left.empty() )
    {
        Range conn_verts;
        rval = context.MBI->get_connectivity( adj_ents_left, conn_verts );MB_CHK_ERR( rval );
        vertices = subtract( vertices, conn_verts );
    }

    rval = context.MBI->delete_entities( vertices );MB_CHK_ERR( rval );

    std::map< std::string, int >::iterator mit;

    for( mit = context.appIdMap.begin(); mit != context.appIdMap.end(); mit++ )
    {
        int pidx = mit->second;

        if( *pid == pidx )
        {
            break;
        }
    }

    context.appIdMap.erase( mit );
    std::map< int, int >::iterator mit1;

    for( mit1 = context.appIdCompMap.begin(); mit1 != context.appIdCompMap.end(); mit1++ )
    {
        int pidx = mit1->second;

        if( *pid == pidx )
        {
            break;
        }
    }

    context.appIdCompMap.erase( mit1 );

    context.unused_pid--;  // we have to go backwards always TODO
    context.appDatas.pop_back();
#ifdef MOAB_HAVE_MPI
    context.pcomms.pop_back();
#endif
    return moab::MB_SUCCESS;
}

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 428 of file iMOAB.cpp.

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

{
    // release any Fortran specific allocations here before we pass it on
    context.appDatas[*pid].is_fortran = false;

    // release all datastructure allocations
    return iMOAB_DeregisterApplication( pid );
}

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 184 of file iMOAB.cpp.

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

Referenced by commgraphtest(), main(), migrate(), and migrate_smart().

{
    context.refCountMB--;

    if( 0 == context.refCountMB )
    {
        delete context.MBI;
    }

    return MB_SUCCESS;
}
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 1140 of file iMOAB.cpp.

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.

Referenced by main().

{
    assert( global_block_ID );      // ensure global block ID argument is not null
    assert( connectivity_length );  // ensure connectivity length argument is not null

    appData& data                     = context.appDatas[*pid];
    std::map< int, int >& matMap      = data.matIndex;
    std::map< int, int >::iterator it = matMap.find( *global_block_ID );

    if( it == matMap.end() )
    {
        return moab::MB_FAILURE;
    }  // error in finding block with id

    int blockIndex          = matMap[*global_block_ID];
    EntityHandle matMeshSet = data.mat_sets[blockIndex];
    std::vector< EntityHandle > elems;

    ErrorCode rval = context.MBI->get_entities_by_handle( matMeshSet, elems );MB_CHK_ERR( rval );

    if( elems.empty() )
    {
        return moab::MB_FAILURE;
    }

    std::vector< EntityHandle > vconnect;
    rval = context.MBI->get_connectivity( &elems[0], elems.size(), vconnect );MB_CHK_ERR( rval );

    if( *connectivity_length != (int)vconnect.size() )
    {
        return moab::MB_FAILURE;
    }  // mismatched sizes

    for( int i = 0; i < *connectivity_length; i++ )
    {
        int inx = data.all_verts.index( vconnect[i] );

        if( -1 == inx )
        {
            return moab::MB_FAILURE;
        }  // error, vertex not in local range

        element_connectivity[i] = inx;
    }

    return moab::MB_SUCCESS;
}
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 1017 of file iMOAB.cpp.

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

Referenced by main().

{
    // local id blocks? they are counted from 0 to number of visible blocks ...
    // will actually return material set tag value for global
    Range& matSets = context.appDatas[*pid].mat_sets;

    if( *block_length != (int)matSets.size() )
    {
        return moab::MB_FAILURE;
    }

    // return material set tag gtags[0 is material set tag
    ErrorCode rval = context.MBI->tag_get_data( context.material_tag, matSets, global_block_IDs );MB_CHK_ERR( rval );

    // populate map with index
    std::map< int, int >& matIdx = context.appDatas[*pid].matIndex;
    for( unsigned i = 0; i < matSets.size(); i++ )
    {
        matIdx[global_block_IDs[i]] = i;
    }

    return moab::MB_SUCCESS;
}
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 1041 of file iMOAB.cpp.

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().

Referenced by main().

{
    assert( global_block_ID );

    std::map< int, int >& matMap      = context.appDatas[*pid].matIndex;
    std::map< int, int >::iterator it = matMap.find( *global_block_ID );

    if( it == matMap.end() )
    {
        return moab::MB_FAILURE;
    }  // error in finding block with id

    int blockIndex          = matMap[*global_block_ID];
    EntityHandle matMeshSet = context.appDatas[*pid].mat_sets[blockIndex];
    Range blo_elems;
    ErrorCode rval = context.MBI->get_entities_by_handle( matMeshSet, blo_elems );

    if( MB_SUCCESS != rval || blo_elems.empty() )
    {
        return moab::MB_FAILURE;
    }

    EntityType type = context.MBI->type_from_handle( blo_elems[0] );

    if( !blo_elems.all_of_type( type ) )
    {
        return moab::MB_FAILURE;
    }  // not all of same  type

    const EntityHandle* conn = NULL;
    int num_verts            = 0;
    rval                     = context.MBI->get_connectivity( blo_elems[0], conn, num_verts );MB_CHK_ERR( rval );

    *vertices_per_element  = num_verts;
    *num_elements_in_block = (int)blo_elems.size();

    return moab::MB_SUCCESS;
}
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 2007 of file iMOAB.cpp.

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.

Referenced by main().

{
    ErrorCode rval;
    // exactly the same code, except tag type check
    std::string tag_names( tag_storage_names );
    // exactly the same code as for int tag :) maybe should check the type of tag too
    std::vector< std::string > tagNames;
    std::vector< Tag > tagHandles;
    std::string separator( ":" );
    split_tag_names( tag_names, separator, tagNames );

    appData& data = context.appDatas[*pid];

    // set it on a subset of entities, based on type and length
    Range* ents_to_get = NULL;

    if( *ent_type == 0 )  // vertices
    {
        ents_to_get = &data.all_verts;
    }
    else if( *ent_type == 1 )
    {
        ents_to_get = &data.primary_elems;
    }
    int nents_to_get = (int)ents_to_get->size();
    int position     = 0;
    for( size_t i = 0; i < tagNames.size(); i++ )
    {
        if( data.tagMap.find( tagNames[i] ) == data.tagMap.end() )
        {
            return moab::MB_FAILURE;
        }  // tag not defined

        Tag tag = data.tagMap[tagNames[i]];

        int tagLength = 0;
        rval          = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval );

        DataType dtype;
        rval = context.MBI->tag_get_data_type( tag, dtype );MB_CHK_ERR( rval );

        if( dtype != MB_TYPE_DOUBLE )
        {
            return moab::MB_FAILURE;
        }

        if( position + nents_to_get * tagLength > *num_tag_storage_length )
            return moab::MB_FAILURE;  // too many entity values to get

        rval = context.MBI->tag_get_data( tag, *ents_to_get, &tag_storage_data[position] );MB_CHK_ERR( rval );
        position = position + nents_to_get * tagLength;
    }

    return moab::MB_SUCCESS;  // no error
}
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 1191 of file iMOAB.cpp.

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().

Referenced by main().

{
    assert( elem_index );           // ensure element index argument is not null
    assert( connectivity_length );  // ensure connectivity length argument is not null

    appData& data = context.appDatas[*pid];
    assert( ( *elem_index >= 0 ) && ( *elem_index < (int)data.primary_elems.size() ) );

    int num_nodes;
    const EntityHandle* conn;

    EntityHandle eh = data.primary_elems[*elem_index];

    ErrorCode rval = context.MBI->get_connectivity( eh, conn, num_nodes );MB_CHK_ERR( rval );

    if( *connectivity_length < num_nodes )
    {
        return moab::MB_FAILURE;
    }  // wrong number of vertices

    for( int i = 0; i < num_nodes; i++ )
    {
        int index = data.all_verts.index( conn[i] );

        if( -1 == index )
        {
            return moab::MB_FAILURE;
        }

        element_connectivity[i] = index;
    }

    *connectivity_length = num_nodes;

    return moab::MB_SUCCESS;
}
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 1281 of file iMOAB.cpp.

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().

Referenced by main().

{
    assert( global_block_ID );        // ensure global block ID argument is not null
    assert( num_elements_in_block );  // ensure number of elements in block argument is not null

    appData& data                = context.appDatas[*pid];
    std::map< int, int >& matMap = data.matIndex;

    std::map< int, int >::iterator it = matMap.find( *global_block_ID );

    if( it == matMap.end() )
    {
        return moab::MB_FAILURE;
    }  // error in finding block with id

    int blockIndex          = matMap[*global_block_ID];
    EntityHandle matMeshSet = data.mat_sets[blockIndex];
    Range elems;
    ErrorCode rval = context.MBI->get_entities_by_handle( matMeshSet, elems );MB_CHK_ERR( rval );

    if( elems.empty() )
    {
        return moab::MB_FAILURE;
    }

    if( *num_elements_in_block != (int)elems.size() )
    {
        return moab::MB_FAILURE;
    }  // bad memory allocation

    rval = context.MBI->tag_get_data( context.globalID_tag, elems, global_element_ID );MB_CHK_ERR( rval );

    // check that elems are among primary_elems in data
    for( int i = 0; i < *num_elements_in_block; i++ )
    {
        local_element_ID[i] = data.primary_elems.index( elems[i] );

        if( -1 == local_element_ID[i] )
        {
            return moab::MB_FAILURE;
        }  // error, not in local primary elements
    }

    return moab::MB_SUCCESS;
}
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 1231 of file iMOAB.cpp.

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, and moab::Range::size().

Referenced by main().

{
    assert( global_block_ID );        // ensure global block ID argument is not null
    assert( num_elements_in_block );  // ensure number of elements in block argument is not null

    std::map< int, int >& matMap = context.appDatas[*pid].matIndex;

    std::map< int, int >::iterator it = matMap.find( *global_block_ID );

    if( it == matMap.end() )
    {
        return moab::MB_FAILURE;
    }  // error in finding block with id

    int blockIndex          = matMap[*global_block_ID];
    EntityHandle matMeshSet = context.appDatas[*pid].mat_sets[blockIndex];
    Range elems;

    ErrorCode rval = context.MBI->get_entities_by_handle( matMeshSet, elems );MB_CHK_ERR( rval );

    if( elems.empty() )
    {
        return moab::MB_FAILURE;
    }

    if( *num_elements_in_block != (int)elems.size() )
    {
        return moab::MB_FAILURE;
    }  // bad memory allocation

    int i = 0;
#ifdef MOAB_HAVE_MPI
    ParallelComm* pco = context.pcomms[*pid];
#endif

    for( Range::iterator vit = elems.begin(); vit != elems.end(); vit++, i++ )
    {
#ifdef MOAB_HAVE_MPI
        rval = pco->get_owner( *vit, element_ownership[i] );MB_CHK_ERR( rval );
#else
        element_ownership[i] = 0; /* owned by 0 */
#endif
    }

    return moab::MB_SUCCESS;
}
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 2284 of file iMOAB.cpp.

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

{
    appData& data = context.appDatas[*pid];
    if( NULL != num_global_verts )
    {
        *num_global_verts = data.num_global_vertices;
    }
    if( NULL != num_global_elems )
    {
        *num_global_elems = data.num_global_elements;
    }

    return moab::MB_SUCCESS;
}
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 1639 of file iMOAB.cpp.

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.

Referenced by main().

{
    ErrorCode rval;
    std::string tag_name( tag_storage_name );

    appData& data = context.appDatas[*pid];

    if( data.tagMap.find( tag_name ) == data.tagMap.end() )
    {
        return moab::MB_FAILURE;
    }  // tag not defined

    Tag tag = data.tagMap[tag_name];

    int tagLength = 0;
    rval          = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval );

    DataType dtype;
    rval = context.MBI->tag_get_data_type( tag, dtype );MB_CHK_ERR( rval );

    if( dtype != MB_TYPE_INTEGER )
    {
        return moab::MB_FAILURE;
    }

    // set it on a subset of entities, based on type and length
    Range* ents_to_get;

    if( *ent_type == 0 )  // vertices
    {
        ents_to_get = &data.all_verts;
    }
    else  // if (*ent_type == 1)
    {
        ents_to_get = &data.primary_elems;
    }

    int nents_to_get = *num_tag_storage_length / tagLength;

    if( nents_to_get > (int)ents_to_get->size() || nents_to_get < 1 )
    {
        return moab::MB_FAILURE;
    }  // to many entities to get, or too little

    // restrict the range; everything is contiguous; or not?
    // Range contig_range ( * ( ents_to_get->begin() ), * ( ents_to_get->begin() + nents_to_get - 1
    // ) );

    rval = context.MBI->tag_get_data( tag, *ents_to_get, tag_storage_data );MB_CHK_ERR( rval );

    return moab::MB_SUCCESS;  // no error
}
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 859 of file iMOAB.cpp.

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.

Referenced by main().

{
    ErrorCode rval;
    appData& data        = context.appDatas[*pid];
    EntityHandle fileSet = data.file_set;

    // this will include ghost elements
    // first clear all data ranges; this can be called after ghosting
    if( num_visible_elements )
    {
        num_visible_elements[2] = static_cast< int >( data.primary_elems.size() );
        // separate ghost and local/owned primary elements
        num_visible_elements[0] = static_cast< int >( data.owned_elems.size() );
        num_visible_elements[1] = static_cast< int >( data.ghost_elems.size() );
    }
    if( num_visible_vertices )
    {
        num_visible_vertices[2] = static_cast< int >( data.all_verts.size() );
        num_visible_vertices[1] = static_cast< int >( data.ghost_vertices.size() );
        // local are those that are not ghosts
        num_visible_vertices[0] = num_visible_vertices[2] - num_visible_vertices[1];
    }

    if( num_visible_blocks )
    {
        rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.material_tag ), 0, 1,
                                                          data.mat_sets, Interface::UNION );MB_CHK_ERR( rval );

        num_visible_blocks[2] = data.mat_sets.size();
        num_visible_blocks[0] = num_visible_blocks[2];
        num_visible_blocks[1] = 0;
    }

    if( num_visible_surfaceBC )
    {
        rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.neumann_tag ), 0, 1,
                                                          data.neu_sets, Interface::UNION );MB_CHK_ERR( rval );

        num_visible_surfaceBC[2] = 0;
        // count how many faces are in each neu set, and how many regions are
        // adjacent to them;
        int numNeuSets = (int)data.neu_sets.size();

        for( int i = 0; i < numNeuSets; i++ )
        {
            Range subents;
            EntityHandle nset = data.neu_sets[i];
            rval              = context.MBI->get_entities_by_dimension( nset, data.dimension - 1, subents );MB_CHK_ERR( rval );

            for( Range::iterator it = subents.begin(); it != subents.end(); ++it )
            {
                EntityHandle subent = *it;
                Range adjPrimaryEnts;
                rval = context.MBI->get_adjacencies( &subent, 1, data.dimension, false, adjPrimaryEnts );MB_CHK_ERR( rval );

                num_visible_surfaceBC[2] += (int)adjPrimaryEnts.size();
            }
        }

        num_visible_surfaceBC[0] = num_visible_surfaceBC[2];
        num_visible_surfaceBC[1] = 0;
    }

    if( num_visible_vertexBC )
    {
        rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.dirichlet_tag ), 0, 1,
                                                          data.diri_sets, Interface::UNION );MB_CHK_ERR( rval );

        num_visible_vertexBC[2] = 0;
        int numDiriSets         = (int)data.diri_sets.size();

        for( int i = 0; i < numDiriSets; i++ )
        {
            Range verts;
            EntityHandle diset = data.diri_sets[i];
            rval               = context.MBI->get_entities_by_dimension( diset, 0, verts );MB_CHK_ERR( rval );

            num_visible_vertexBC[2] += (int)verts.size();
        }

        num_visible_vertexBC[0] = num_visible_vertexBC[2];
        num_visible_vertexBC[1] = 0;
    }

    return moab::MB_SUCCESS;
}
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 2152 of file iMOAB.cpp.

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().

Referenced by main().

{
    ErrorCode rval;

    // one neighbor for each subentity of dimension-1
    MeshTopoUtil mtu( context.MBI );
    appData& data   = context.appDatas[*pid];
    EntityHandle eh = data.primary_elems[*local_index];
    Range adjs;
    rval = mtu.get_bridge_adjacencies( eh, data.dimension - 1, data.dimension, adjs );MB_CHK_ERR( rval );

    if( *num_adjacent_elements < (int)adjs.size() )
    {
        return moab::MB_FAILURE;
    }  // not dimensioned correctly

    *num_adjacent_elements = (int)adjs.size();

    for( int i = 0; i < *num_adjacent_elements; i++ )
    {
        adjacent_element_IDs[i] = data.primary_elems.index( adjs[i] );
    }

    return moab::MB_SUCCESS;
}
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 1331 of file iMOAB.cpp.

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::side_number(), moab::Interface::side_number(), moab::Range::size(), and moab::Interface::tag_get_data().

Referenced by main().

{
    // we have to fill bc data for neumann sets;/
    ErrorCode rval;

    // it was counted above, in GetMeshInfo
    appData& data  = context.appDatas[*pid];
    int numNeuSets = (int)data.neu_sets.size();

    int index = 0;  // index [0, surface_BC_length) for the arrays returned

    for( int i = 0; i < numNeuSets; i++ )
    {
        Range subents;
        EntityHandle nset = data.neu_sets[i];
        rval              = context.MBI->get_entities_by_dimension( nset, data.dimension - 1, subents );MB_CHK_ERR( rval );

        int neuVal;
        rval = context.MBI->tag_get_data( context.neumann_tag, &nset, 1, &neuVal );MB_CHK_ERR( rval );

        for( Range::iterator it = subents.begin(); it != subents.end(); ++it )
        {
            EntityHandle subent = *it;
            Range adjPrimaryEnts;
            rval = context.MBI->get_adjacencies( &subent, 1, data.dimension, false, adjPrimaryEnts );MB_CHK_ERR( rval );

            // get global id of the primary ents, and side number of the quad/subentity
            // this is moab ordering
            for( Range::iterator pit = adjPrimaryEnts.begin(); pit != adjPrimaryEnts.end(); pit++ )
            {
                EntityHandle primaryEnt = *pit;
                // get global id
                /*int globalID;
                rval = context.MBI->tag_get_data(gtags[3], &primaryEnt, 1, &globalID);
                if (MB_SUCCESS!=rval)
                  return moab::MB_FAILURE;
                global_element_ID[index] = globalID;*/

                // get local element id
                local_element_ID[index] = data.primary_elems.index( primaryEnt );

                if( -1 == local_element_ID[index] )
                {
                    return moab::MB_FAILURE;
                }  // did not find the element locally

                int side_number, sense, offset;
                rval = context.MBI->side_number( primaryEnt, subent, side_number, sense, offset );MB_CHK_ERR( rval );

                reference_surface_ID[index]     = side_number + 1;  // moab is from 0 to 5, it needs 1 to 6
                boundary_condition_value[index] = neuVal;
                index++;
            }
        }
    }

    if( index != *surface_BC_length )
    {
        return moab::MB_FAILURE;
    }  // error in array allocations

    return moab::MB_SUCCESS;
}
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 1399 of file iMOAB.cpp.

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().

Referenced by main().

{
    ErrorCode rval;

    // it was counted above, in GetMeshInfo
    appData& data   = context.appDatas[*pid];
    int numDiriSets = (int)data.diri_sets.size();
    int index       = 0;  // index [0, *vertex_BC_length) for the arrays returned

    for( int i = 0; i < numDiriSets; i++ )
    {
        Range verts;
        EntityHandle diset = data.diri_sets[i];
        rval               = context.MBI->get_entities_by_dimension( diset, 0, verts );MB_CHK_ERR( rval );

        int diriVal;
        rval = context.MBI->tag_get_data( context.dirichlet_tag, &diset, 1, &diriVal );MB_CHK_ERR( rval );

        for( Range::iterator vit = verts.begin(); vit != verts.end(); ++vit )
        {
            EntityHandle vt = *vit;
            /*int vgid;
            rval = context.MBI->tag_get_data(gtags[3], &vt, 1, &vgid);
            if (MB_SUCCESS!=rval)
              return moab::MB_FAILURE;
            global_vertext_ID[index] = vgid;*/
            local_vertex_ID[index] = data.all_verts.index( vt );

            if( -1 == local_vertex_ID[index] )
            {
                return moab::MB_FAILURE;
            }  // vertex was not found

            boundary_condition_value[index] = diriVal;
            index++;
        }
    }

    if( *vertex_BC_length != index )
    {
        return moab::MB_FAILURE;
    }  // array allocation issue

    return moab::MB_SUCCESS;
}
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 951 of file iMOAB.cpp.

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

Referenced by main().

{
    IMOAB_CHECKPOINTER( vertices_length, 2 );
    IMOAB_CHECKPOINTER( global_vertex_ID, 3 );

    const Range& verts = context.appDatas[*pid].all_verts;
    // check for problems with array length
    IMOAB_ASSERT( *vertices_length == static_cast< int >( verts.size() ), "Invalid vertices length provided" );

    // global id tag is context.globalID_tag
    return context.MBI->tag_get_data( context.globalID_tag, verts, global_vertex_ID );
}
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 964 of file iMOAB.cpp.

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

Referenced by main().

{
    assert( vertices_length && *vertices_length );
    assert( visible_global_rank_ID );

    Range& verts = context.appDatas[*pid].all_verts;
    int i        = 0;
#ifdef MOAB_HAVE_MPI
    ParallelComm* pco = context.pcomms[*pid];

    for( Range::iterator vit = verts.begin(); vit != verts.end(); vit++, i++ )
    {
        ErrorCode rval = pco->get_owner( *vit, visible_global_rank_ID[i] );MB_CHK_ERR( rval );
    }

    if( i != *vertices_length )
    {
        return moab::MB_FAILURE;
    }  // warning array allocation problem

#else

    /* everything owned by proc 0 */
    if( (int)verts.size() != *vertices_length )
    {
        return moab::MB_FAILURE;
    }  // warning array allocation problem

    for( Range::iterator vit = verts.begin(); vit != verts.end(); vit++, i++ )
    {
        visible_global_rank_ID[i] = 0;
    }  // all vertices are owned by processor 0, as this is serial run

#endif

    return moab::MB_SUCCESS;
}
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.

Todo:
: Can we also return the index for each block ID?
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 1083 of file iMOAB.cpp.

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, appData::primary_elems, and moab::Interface::tag_get_data().

Referenced by main().

{
    appData& data = context.appDatas[*pid];
#ifdef MOAB_HAVE_MPI
    ParallelComm* pco = context.pcomms[*pid];
#endif

    ErrorCode rval = context.MBI->tag_get_data( context.globalID_tag, data.primary_elems, element_global_IDs );MB_CHK_ERR( rval );

    int i = 0;

    for( Range::iterator eit = data.primary_elems.begin(); eit != data.primary_elems.end(); ++eit, ++i )
    {
#ifdef MOAB_HAVE_MPI
        rval = pco->get_owner( *eit, ranks[i] );MB_CHK_ERR( rval );

#else
        /* everything owned by task 0 */
        ranks[i]             = 0;
#endif
    }

    for( Range::iterator mit = data.mat_sets.begin(); mit != data.mat_sets.end(); ++mit )
    {
        EntityHandle matMeshSet = *mit;
        Range elems;
        rval = context.MBI->get_entities_by_handle( matMeshSet, elems );MB_CHK_ERR( rval );

        int valMatTag;
        rval = context.MBI->tag_get_data( context.material_tag, &matMeshSet, 1, &valMatTag );MB_CHK_ERR( rval );

        for( Range::iterator eit = elems.begin(); eit != elems.end(); ++eit )
        {
            EntityHandle eh = *eit;
            int index       = data.primary_elems.index( eh );

            if( -1 == index )
            {
                return moab::MB_FAILURE;
            }

            if( -1 >= *num_visible_elements )
            {
                return moab::MB_FAILURE;
            }

            block_IDs[index] = valMatTag;
        }
    }

    return moab::MB_SUCCESS;
}
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 1002 of file iMOAB.cpp.

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

Referenced by main().

{
    Range& verts = context.appDatas[*pid].all_verts;

    // interleaved coordinates, so that means deep copy anyway
    if( *coords_length != 3 * (int)verts.size() )
    {
        return moab::MB_FAILURE;
    }

    ErrorCode rval = context.MBI->get_coords( verts, coordinates );MB_CHK_ERR( rval );

    return moab::MB_SUCCESS;
}
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 133 of file iMOAB.cpp.

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, MPI_COMM_WORLD, GlobalContext::MPI_initialized, NEUMANN_SET_TAG_NAME, GlobalContext::neumann_tag, GlobalContext::refCountMB, moab::Interface::tag_get_handle(), and GlobalContext::worldprocs.

Referenced by commgraphtest(), iMOAB_InitializeFortran(), main(), migrate(), and migrate_smart().

{
    if( argc ) IMOAB_CHECKPOINTER( argv, 1 );

    context.iArgc = argc;
    context.iArgv = argv;  // shallow copy

    if( 0 == context.refCountMB )
    {
        context.MBI = new( std::nothrow ) moab::Core;
        // retrieve the default tags
        const char* const shared_set_tag_names[] = { MATERIAL_SET_TAG_NAME, NEUMANN_SET_TAG_NAME,
                                                     DIRICHLET_SET_TAG_NAME, GLOBAL_ID_TAG_NAME };
        // blocks, visible surfaceBC(neumann), vertexBC (Dirichlet), global id, parallel partition
        Tag gtags[4];

        for( int i = 0; i < 4; i++ )
        {

            ErrorCode rval =
                context.MBI->tag_get_handle( shared_set_tag_names[i], 1, MB_TYPE_INTEGER, gtags[i], MB_TAG_ANY );MB_CHK_ERR( rval );
        }

        context.material_tag  = gtags[0];
        context.neumann_tag   = gtags[1];
        context.dirichlet_tag = gtags[2];
        context.globalID_tag  = gtags[3];
    }

    context.MPI_initialized = false;
#ifdef MOAB_HAVE_MPI
    int flagInit;
    MPI_Initialized( &flagInit );

    if( flagInit && !context.MPI_initialized )
    {
        MPI_Comm_size( MPI_COMM_WORLD, &context.worldprocs );
        MPI_Comm_rank( MPI_COMM_WORLD, &context.globalrank );
        context.MPI_initialized = true;
    }
#endif

    context.refCountMB++;
    return moab::MB_SUCCESS;
}

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 179 of file iMOAB.cpp.

References iMOAB_Initialize().

{
    return iMOAB_Initialize( 0, 0 );
}
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 603 of file iMOAB.cpp.

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

Referenced by commgraphtest(), main(), migrate(), migrate_smart(), and setup_component_coupler_meshes().

{
    IMOAB_CHECKPOINTER( filename, 2 );
    IMOAB_ASSERT( strlen( filename ), "Invalid filename length." );

    // make sure we use the file set and pcomm associated with the *pid
    std::ostringstream newopts;
    if( read_options ) newopts << read_options;

#ifdef MOAB_HAVE_MPI

    if( context.MPI_initialized )
    {
        if( context.worldprocs > 1 )
        {
            std::string opts( ( read_options ? read_options : "" ) );
            std::string pcid( "PARALLEL_COMM=" );
            std::size_t found = opts.find( pcid );

            if( found != std::string::npos )
            {
                std::cerr << " cannot specify PARALLEL_COMM option, it is implicit \n";
                return moab::MB_FAILURE;
            }

            // in serial, apply PARALLEL_COMM option only for h5m files; it does not work for .g
            // files (used in test_remapping)
            std::string filen( filename );
            std::string::size_type idx = filen.rfind( '.' );

            if( idx != std::string::npos )
            {
                std::string extension = filen.substr( idx + 1 );
                if( (extension == std::string( "h5m" )) || (extension == std::string( "nc" ))) newopts << ";;PARALLEL_COMM=" << *pid;
            }

            if( *num_ghost_layers >= 1 )
            {
                // if we want ghosts, we will want additional entities, the last .1
                // because the addl ents can be edges, faces that are part of the neumann sets
                std::string pcid2( "PARALLEL_GHOSTS=" );
                std::size_t found2 = opts.find( pcid2 );

                if( found2 != std::string::npos )
                {
                    std::cout << " PARALLEL_GHOSTS option is already specified, ignore passed "
                                 "number of layers \n";
                }
                else
                {
                    // dimension of primary entities is 3 here, but it could be 2 for climate
                    // meshes; we would need to pass PARALLEL_GHOSTS explicitly for 2d meshes, for
                    // example:  ";PARALLEL_GHOSTS=2.0.1"
                    newopts << ";PARALLEL_GHOSTS=3.0." << *num_ghost_layers << ".3";
                }
            }
        }
    }
#else
    IMOAB_ASSERT( *num_ghost_layers == 0, "Cannot provide ghost layers in serial." );
#endif

    // Now let us actually load the MOAB file with the appropriate read options
    ErrorCode rval = context.MBI->load_file( filename, &context.appDatas[*pid].file_set, newopts.str().c_str() );MB_CHK_ERR( rval );

#ifdef VERBOSE
    // some debugging stuff
    std::ostringstream outfile;
#ifdef MOAB_HAVE_MPI
    int rank   = context.pcomms[*pid]->rank();
    int nprocs = context.pcomms[*pid]->size();
    outfile << "TaskMesh_n" << nprocs << "." << rank << ".h5m";
#else
    outfile << "TaskMesh_n1.0.h5m";
#endif
    // the mesh contains ghosts too, but they are not part of mat/neumann set
    // write in serial the file, to see what tags are missing
    rval = context.MBI->write_file( outfile.str().c_str() );MB_CHK_ERR( rval );  // everything on current task, written in serial
#endif

    // Update mesh information
    return iMOAB_UpdateMeshInfo( pid );
}
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 459 of file iMOAB.cpp.

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.

Referenced by fdriver(), and main().

{
    IMOAB_CHECKPOINTER( filename, 1 );
    IMOAB_ASSERT( strlen( filename ), "Invalid filename length." );

#ifdef MOAB_HAVE_HDF5
    std::string filen( filename );

    int edges   = 0;
    int faces   = 0;
    int regions = 0;
    if( num_global_vertices ) *num_global_vertices = 0;
    if( num_global_elements ) *num_global_elements = 0;
    if( num_dimension ) *num_dimension = 0;
    if( num_parts ) *num_parts = 0;

    mhdf_FileHandle file;
    mhdf_Status status;
    unsigned long max_id;
    struct mhdf_FileDesc* data;

    file = mhdf_openFile( filen.c_str(), 0, &max_id, -1, &status );

    if( mhdf_isError( &status ) )
    {
        fprintf( stderr, "%s: %s\n", filename, mhdf_message( &status ) );
        return moab::MB_FAILURE;
    }

    data = mhdf_getFileSummary( file, H5T_NATIVE_ULONG, &status,
                                1 );  // will use extra set info; will get parallel partition tag info too!

    if( mhdf_isError( &status ) )
    {
        fprintf( stderr, "%s: %s\n", filename, mhdf_message( &status ) );
        return moab::MB_FAILURE;
    }

    if( num_dimension ) *num_dimension = data->nodes.vals_per_ent;
    if( num_global_vertices ) *num_global_vertices = (int)data->nodes.count;

    for( int i = 0; i < data->num_elem_desc; i++ )
    {
        struct mhdf_ElemDesc* el_desc = &( data->elems[i] );
        struct mhdf_EntDesc* ent_d    = &( el_desc->desc );

        if( 0 == strcmp( el_desc->type, mhdf_EDGE_TYPE_NAME ) )
        {
            edges += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mhdf_TRI_TYPE_NAME ) )
        {
            faces += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mhdf_QUAD_TYPE_NAME ) )
        {
            faces += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mhdf_POLYGON_TYPE_NAME ) )
        {
            faces += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mhdf_TET_TYPE_NAME ) )
        {
            regions += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mhdf_PYRAMID_TYPE_NAME ) )
        {
            regions += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mhdf_PRISM_TYPE_NAME ) )
        {
            regions += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mdhf_KNIFE_TYPE_NAME ) )
        {
            regions += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mdhf_HEX_TYPE_NAME ) )
        {
            regions += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mhdf_POLYHEDRON_TYPE_NAME ) )
        {
            regions += ent_d->count;
        }

        if( 0 == strcmp( el_desc->type, mhdf_SEPTAHEDRON_TYPE_NAME ) )
        {
            regions += ent_d->count;
        }
    }

    if( num_parts ) *num_parts = data->numEntSets[0];

    // is this required?
    if( edges > 0 )
    {
        if( num_dimension ) *num_dimension = 1;  // I don't think it will ever return 1
        if( num_global_elements ) *num_global_elements = edges;
    }

    if( faces > 0 )
    {
        if( num_dimension ) *num_dimension = 2;
        if( num_global_elements ) *num_global_elements = faces;
    }

    if( regions > 0 )
    {
        if( num_dimension ) *num_dimension = 3;
        if( num_global_elements ) *num_global_elements = regions;
    }

    mhdf_closeFile( file, &status );

    free( data );

#else
    std::cout << filename
              << ": Please reconfigure with HDF5. Cannot retrieve header information for file "
                 "formats other than a h5m file.\n";
    if( num_global_vertices ) *num_global_vertices = 0;
    if( num_global_elements ) *num_global_elements = 0;
    if( num_dimension ) *num_dimension = 0;
    if( num_parts ) *num_parts = 0;
#endif

    return moab::MB_SUCCESS;
}
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 2112 of file iMOAB.cpp.

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

Referenced by main().

{
#ifdef MOAB_HAVE_MPI
    appData& data = context.appDatas[*pid];
    Range ent_exchange;

    if( *tag_index < 0 || *tag_index >= (int)data.tagList.size() )
    {
        return moab::MB_FAILURE;
    }  // error in tag index

    Tag tagh = data.tagList[*tag_index];

    if( *ent_type == 0 )
    {
        ent_exchange = data.all_verts;
    }
    else if( *ent_type == 1 )
    {
        ent_exchange = data.primary_elems;
    }
    else
    {
        return moab::MB_FAILURE;
    }  // unexpected type

    ParallelComm* pco = context.pcomms[*pid];
    // we could do different MPI_Op; do not bother now, we will call from fortran
    ErrorCode rval = pco->reduce_tags( tagh, MPI_MAX, ent_exchange );MB_CHK_ERR( rval );

#else
    /* do nothing if serial */
    // just silence the warning
    // do not call sync tags in serial!
    int k = *pid + *tag_index + *ent_type;
    k++;  // just do junk, to avoid complaints
#endif
    return moab::MB_SUCCESS;
}
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.

Definition at line 196 of file iMOAB.cpp.

References GlobalContext::appDatas, GlobalContext::appIdCompMap, GlobalContext::appIdMap, context, moab::Interface::create_meshset(), ErrorCode, appData::file_set, moab::ParallelComm::get_id(), appData::global_id, IMOAB_CHECKPOINTER, appData::is_fortran, MB_CHK_ERR, MB_SUCCESS, GlobalContext::MBI, MESHSET_SET, appData::name, appData::point_cloud, and GlobalContext::unused_pid.

Referenced by commgraphtest(), iMOAB_RegisterApplicationFortran(), main(), migrate(), and migrate_smart().

{
    IMOAB_CHECKPOINTER( app_name, 1 );
#ifdef MOAB_HAVE_MPI
    IMOAB_CHECKPOINTER( comm, 2 );
    IMOAB_CHECKPOINTER( compid, 3 );
#else
    IMOAB_CHECKPOINTER( compid, 2 );
#endif

    // will create a parallel comm for this application too, so there will be a
    // mapping from *pid to file set and to parallel comm instances
    std::string name( app_name );

    if( context.appIdMap.find( name ) != context.appIdMap.end() )
    {
        std::cout << " application " << name << " already registered \n";
        return moab::MB_FAILURE;
    }

    *pid                   = context.unused_pid++;
    context.appIdMap[name] = *pid;
    int rankHere           = 0;
#ifdef MOAB_HAVE_MPI
    MPI_Comm_rank( *comm, &rankHere );
#endif
    if( !rankHere ) std::cout << " application " << name << " with ID = " << *pid << " is registered now \n";
    if( *compid <= 0 )
    {
        std::cout << " convention for external application is to have its id positive \n";
        return moab::MB_FAILURE;
    }

    if( context.appIdCompMap.find( *compid ) != context.appIdCompMap.end() )
    {
        std::cout << " external application with comp id " << *compid << " is already registered\n";
        return moab::MB_FAILURE;
    }

    context.appIdCompMap[*compid] = *pid;

    // now create ParallelComm and a file set for this application
#ifdef MOAB_HAVE_MPI
    if( *comm )
    {
        ParallelComm* pco = new ParallelComm( context.MBI, *comm );

#ifndef NDEBUG
        int index = pco->get_id();  // it could be useful to get app id from pcomm instance ...
        assert( index == *pid );
        // here, we assert the the pid is the same as the id of the ParallelComm instance
        // useful for writing in parallel
#endif
        context.pcomms.push_back( pco );
    }
    else
    {
        context.pcomms.push_back( 0 );
    }
#endif

    // create now the file set that will be used for loading the model in
    EntityHandle file_set;
    ErrorCode rval = context.MBI->create_meshset( MESHSET_SET, file_set );MB_CHK_ERR( rval );

    appData app_data;
    app_data.file_set  = file_set;
    app_data.global_id = *compid;  // will be used mostly for par comm graph
    app_data.name      = name;     // save the name of application

#ifdef MOAB_HAVE_TEMPESTREMAP
    app_data.tempestData.remapper = NULL;  // Only allocate as needed
#endif

    app_data.point_cloud = false;
    app_data.is_fortran  = false;

    context.appDatas.push_back(
        app_data );  // it will correspond to app_FileSets[*pid] will be the file set of interest
    return moab::MB_SUCCESS;
}
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.

Definition at line 283 of file iMOAB.cpp.

References GlobalContext::appDatas, context, ErrCode, IMOAB_CHECKPOINTER, and iMOAB_RegisterApplication().

{
    IMOAB_CHECKPOINTER( app_name, 1 );
#ifdef MOAB_HAVE_MPI
    IMOAB_CHECKPOINTER( comm, 2 );
    IMOAB_CHECKPOINTER( compid, 3 );
#else
    IMOAB_CHECKPOINTER( compid, 2 );
#endif

    ErrCode err;
    assert( app_name != nullptr );
    std::string name( app_name );

#ifdef MOAB_HAVE_MPI
    MPI_Comm ccomm;
    if( comm )
    {
        // convert from Fortran communicator to a C communicator
        // see transfer of handles
        // http://www.mpi-forum.org/docs/mpi-2.2/mpi22-report/node361.htm
        ccomm = MPI_Comm_f2c( (MPI_Fint)*comm );
    }
#endif

    // now call C style registration function:
    err = iMOAB_RegisterApplication( app_name,
#ifdef MOAB_HAVE_MPI
                                     &ccomm,
#endif
                                     compid, pid );

    // Now that we have created the application context, store that
    // the application being registered is from a Fortran context
    context.appDatas[*pid].is_fortran = true;

    return err;
}
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 1696 of file iMOAB.cpp.

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.

Referenced by main().

{
    ErrorCode rval;
    std::string tag_names( tag_storage_names );
    // exactly the same code as for int tag :) maybe should check the type of tag too
    std::vector< std::string > tagNames;
    std::vector< Tag > tagHandles;
    std::string separator( ":" );
    split_tag_names( tag_names, separator, tagNames );

    appData& data      = context.appDatas[*pid];
    Range* ents_to_set = NULL;

    if( *ent_type == 0 )  // vertices
    {
        ents_to_set = &data.all_verts;
    }
    else if( *ent_type == 1 )
    {
        ents_to_set = &data.primary_elems;
    }

    int nents_to_be_set = (int)( *ents_to_set ).size();
    int position        = 0;

    for( size_t i = 0; i < tagNames.size(); i++ )
    {
        if( data.tagMap.find( tagNames[i] ) == data.tagMap.end() )
        {
            return moab::MB_FAILURE;
        }  // some tag not defined yet in the app

        Tag tag = data.tagMap[tagNames[i]];

        int tagLength = 0;
        rval          = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval );

        DataType dtype;
        rval = context.MBI->tag_get_data_type( tag, dtype );MB_CHK_ERR( rval );

        if( dtype != MB_TYPE_DOUBLE )
        {
            return moab::MB_FAILURE;
        }

        // set it on the subset of entities, based on type and length
        if( position + tagLength * nents_to_be_set > *num_tag_storage_length )
            return moab::MB_FAILURE;  // too many entity values to be set

        rval = context.MBI->tag_set_data( tag, *ents_to_set, &tag_storage_data[position] );MB_CHK_ERR( rval );
        // increment position to next tag
        position = position + tagLength * nents_to_be_set;
    }
    return moab::MB_SUCCESS;  // no error
}
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 1756 of file iMOAB.cpp.

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, 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, and moab::TupleList::vr_rd.

Referenced by main().

{
    ErrorCode rval;
    std::string tag_names( tag_storage_names );
    // exactly the same code as for int tag :) maybe should check the type of tag too
    std::vector< std::string > tagNames;
    std::vector< Tag > tagHandles;
    std::string separator( ":" );
    split_tag_names( tag_names, separator, tagNames );

    appData& data      = context.appDatas[*pid];
    Range* ents_to_set = NULL;

    if( *ent_type == 0 )  // vertices
    {
        ents_to_set = &data.all_verts;
    }
    else if( *ent_type == 1 )
    {
        ents_to_set = &data.primary_elems;
    }

    int nents_to_be_set = (int)( *ents_to_set ).size();

    Tag gidTag = context.MBI->globalId_tag();
    std::vector< int > gids;
    gids.resize( nents_to_be_set );
    rval = context.MBI->tag_get_data( gidTag, *ents_to_set, &gids[0] );MB_CHK_ERR( rval );

    // so we will need to set the tags according to the global id passed;
    // so the order in tag_storage_data is the same as the order in globalIds, but the order
    // in local range is gids
    std::map< int, EntityHandle > eh_by_gid;
    int i = 0;
    for( Range::iterator it = ents_to_set->begin(); it != ents_to_set->end(); ++it, ++i )
    {
        eh_by_gid[gids[i]] = *it;
    }

    std::vector< int > tagLengths( tagNames.size() );
    std::vector< Tag > tagList;
    size_t total_tag_len = 0;
    for( size_t i = 0; i < tagNames.size(); i++ )
    {
        if( data.tagMap.find( tagNames[i] ) == data.tagMap.end() )
        {
            MB_SET_ERR( moab::MB_FAILURE, "tag missing" );
        }  // some tag not defined yet in the app

        Tag tag = data.tagMap[tagNames[i]];
        tagList.push_back( tag );

        int tagLength = 0;
        rval          = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval );

        total_tag_len += tagLength;
        tagLengths[i] = tagLength;
        DataType dtype;
        rval = context.MBI->tag_get_data_type( tag, dtype );MB_CHK_ERR( rval );

        if( dtype != MB_TYPE_DOUBLE )
        {
            MB_SET_ERR( moab::MB_FAILURE, "tag not double type" );
        }
    }
    bool serial = true;
#ifdef MOAB_HAVE_MPI
    ParallelComm* pco = context.pcomms[*pid];
    unsigned num_procs     = pco->size();
    if( num_procs > 1 ) serial = false;
#endif

    if( serial )
    {
        assert( total_tag_len * nents_to_be_set - *num_tag_storage_length == 0 );
        // tags are unrolled, we loop over global ids first, then careful about tags
        for( int i = 0; i < nents_to_be_set; i++ )
        {
            int gid         = globalIds[i];
            EntityHandle eh = eh_by_gid[gid];
            // now loop over tags
            int indexInTagValues = 0;  //
            for( size_t j = 0; j < tagList.size(); j++ )
            {
                indexInTagValues += i * tagLengths[j];
                rval = context.MBI->tag_set_data( tagList[j], &eh, 1, &tag_storage_data[indexInTagValues] );MB_CHK_ERR( rval );
                // advance the pointer/index
                indexInTagValues += ( nents_to_be_set - i ) * tagLengths[j];  // at the end of tag data
            }
        }
    }
#ifdef MOAB_HAVE_MPI
    else  // it can be not serial only if pco->size() > 1, parallel
    {
        // in this case, we have to use 2 crystal routers, to send data to the processor that needs it
        // we will create first a tuple to rendevous points, then from there send to the processor that requested it
        // it is a 2-hop global gather scatter
        int nbLocalVals = *num_tag_storage_length / ( (int)tagNames.size() );
        assert( nbLocalVals * tagNames.size() - *num_tag_storage_length == 0 );
        TupleList TLsend;
        TLsend.initialize( 2, 0, 0, total_tag_len, nbLocalVals );  //  to proc, marker(gid), total_tag_len doubles
        TLsend.enableWriteAccess();
        // the processor id that processes global_id is global_id / num_ents_per_proc

        int indexInRealLocal = 0;
        for( int i = 0; i < nbLocalVals; i++ )
        {
            // to proc, marker, element local index, index in el
            int marker              = globalIds[i];
            int to_proc             = marker % num_procs;
            int n                   = TLsend.get_n();
            TLsend.vi_wr[2 * n]     = to_proc;  // send to processor
            TLsend.vi_wr[2 * n + 1] = marker;
            int indexInTagValues = 0;
            // tag data collect by number of tags
            for( size_t j = 0; j < tagList.size(); j++ )
            {
                indexInTagValues += i * tagLengths[j];
                for( int k = 0; k < tagLengths[j]; k++ )
                {
                    TLsend.vr_wr[indexInRealLocal++] = tag_storage_data[indexInTagValues + k];
                }
                indexInTagValues += ( nbLocalVals - i ) * tagLengths[j];
            }
            TLsend.inc_n();
        }

        assert( nbLocalVals * total_tag_len - indexInRealLocal == 0 );
        // send now requests, basically inform the rendez-vous point who needs a particular global id
        // send the data to the other processors:
        ( pco->proc_config().crystal_router() )->gs_transfer( 1, TLsend, 0 );
        TupleList TLreq;
        TLreq.initialize( 2, 0, 0, 0, nents_to_be_set );
        TLreq.enableWriteAccess();
        for( int i = 0; i < nents_to_be_set; i++ )
        {
            // to proc, marker
            int marker             = gids[i];
            int to_proc            = marker % num_procs;
            int n                  = TLreq.get_n();
            TLreq.vi_wr[2 * n]     = to_proc;  // send to processor
            TLreq.vi_wr[2 * n + 1] = marker;
            // tag data collect by number of tags
            TLreq.inc_n();
        }
        ( pco->proc_config().crystal_router() )->gs_transfer( 1, TLreq, 0 );

        // we know now that process TLreq.vi_wr[2 * n] needs tags for gid TLreq.vi_wr[2 * n + 1]
        // we should first order by global id, and then build the new TL with send to proc, global id and
        // tags for it
        // sort by global ids the tuple lists
        moab::TupleList::buffer sort_buffer;
        sort_buffer.buffer_init( TLreq.get_n() );
        TLreq.sort( 1, &sort_buffer );
        sort_buffer.reset();
        sort_buffer.buffer_init( TLsend.get_n() );
        TLsend.sort( 1, &sort_buffer );
        sort_buffer.reset();
        // now send the tag values to the proc that requested it
        // in theory, for a full  partition, TLreq  and TLsend should have the same size, and
        // each dof should have exactly one target proc. Is that true or not in general ?
        // how do we plan to use this? Is it better to store the comm graph for future

        // start copy from comm graph settle
        TupleList TLBack;
        TLBack.initialize( 3, 0, 0, total_tag_len, 0 );  // to proc, marker, tag from proc , tag values
        TLBack.enableWriteAccess();

        int n1 = TLreq.get_n();
        int n2 = TLsend.get_n();

        int indexInTLreq  = 0;
        int indexInTLsend = 0;  // advance both, according to the marker
        if( n1 > 0 && n2 > 0 )
        {

            while( indexInTLreq < n1 && indexInTLsend < n2 )  // if any is over, we are done
            {
                int currentValue1 = TLreq.vi_rd[2 * indexInTLreq + 1];
                int currentValue2 = TLsend.vi_rd[2 * indexInTLsend + 1];
                if( currentValue1 < currentValue2 )
                {
                    // we have a big problem; basically, we are saying that
                    // dof currentValue is on one model and not on the other
                    // std::cout << " currentValue1:" << currentValue1 << " missing in comp2" << "\n";
                    indexInTLreq++;
                    continue;
                }
                if( currentValue1 > currentValue2 )
                {
                    // std::cout << " currentValue2:" << currentValue2 << " missing in comp1" << "\n";
                    indexInTLsend++;
                    continue;
                }
                int size1 = 1;
                int size2 = 1;
                while( indexInTLreq + size1 < n1 && currentValue1 == TLreq.vi_rd[2 * ( indexInTLreq + size1 ) + 1] )
                    size1++;
                while( indexInTLsend + size2 < n2 && currentValue2 == TLsend.vi_rd[2 * ( indexInTLsend + size2 ) + 1] )
                    size2++;
                // must be found in both lists, find the start and end indices
                for( int i1 = 0; i1 < size1; i1++ )
                {
                    for( int i2 = 0; i2 < size2; i2++ )
                    {
                        // send the info back to components
                        int n = TLBack.get_n();
                        TLBack.reserve();
                        TLBack.vi_wr[3 * n] = TLreq.vi_rd[2 * ( indexInTLreq + i1 )];  // send back to the proc marker
                                                                                       // came from, info from comp2
                        TLBack.vi_wr[3 * n + 1] = currentValue1;  // initial value (resend, just for verif ?)
                        TLBack.vi_wr[3 * n + 2] = TLsend.vi_rd[2 * ( indexInTLsend + i2 )];  // from proc on comp2
                        // also fill tag values
                        for( size_t k = 0; k < total_tag_len; k++ )
                        {
                            TLBack.vr_rd[total_tag_len * n + k] =
                                TLsend.vr_rd[total_tag_len * indexInTLsend + k];  // deep copy of tag values
                        }
                    }
                }
                indexInTLreq += size1;
                indexInTLsend += size2;
            }
        }
        ( pco->proc_config().crystal_router() )->gs_transfer( 1, TLBack, 0 );
        // end copy from comm graph
        // after we are done sending, we need to set those tag values, in a reverse process compared to send
        n1             = TLBack.get_n();
        double* ptrVal = &TLBack.vr_rd[0];  //
        for( int i = 0; i < n1; i++ )
        {
            int gid         = TLBack.vi_rd[3 * i + 1];  // marker
            EntityHandle eh = eh_by_gid[gid];
            // now loop over tags

            for( size_t j = 0; j < tagList.size(); j++ )
            {
                rval = context.MBI->tag_set_data( tagList[j], &eh, 1, (void*)ptrVal );MB_CHK_ERR( rval );
                // advance the pointer/index
                ptrVal += tagLengths[j];  // at the end of tag data per call
            }
        }
    }
#endif
    return MB_SUCCESS;
}
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 2276 of file iMOAB.cpp.

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

Referenced by main().

{
    appData& data            = context.appDatas[*pid];
    data.num_global_vertices = *num_global_verts;
    data.num_global_elements = *num_global_elems;
    return moab::MB_SUCCESS;
}
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 1582 of file iMOAB.cpp.

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.

Referenced by main().

{
    std::string tag_name( tag_storage_name );

    // Get the application data
    appData& data = context.appDatas[*pid];

    if( data.tagMap.find( tag_name ) == data.tagMap.end() )
    {
        return moab::MB_FAILURE;
    }  // tag not defined

    Tag tag = data.tagMap[tag_name];

    int tagLength  = 0;
    ErrorCode rval = context.MBI->tag_get_length( tag, tagLength );MB_CHK_ERR( rval );

    DataType dtype;
    rval = context.MBI->tag_get_data_type( tag, dtype );

    if( MB_SUCCESS != rval || dtype != MB_TYPE_INTEGER )
    {
        return moab::MB_FAILURE;
    }

    // set it on a subset of entities, based on type and length
    Range* ents_to_set;

    if( *ent_type == 0 )  // vertices
    {
        ents_to_set = &data.all_verts;
    }
    else  // if (*ent_type == 1) // *ent_type can be 0 (vertices) or 1 (elements)
    {
        ents_to_set = &data.primary_elems;
    }

    int nents_to_be_set = *num_tag_storage_length / tagLength;

    if( nents_to_be_set > (int)ents_to_set->size() || nents_to_be_set < 1 )
    {
        return moab::MB_FAILURE;
    }  // to many entities to be set or too few

    // restrict the range; everything is contiguous; or not?

    // Range contig_range ( * ( ents_to_set->begin() ), * ( ents_to_set->begin() + nents_to_be_set -
    // 1 ) );
    rval = context.MBI->tag_set_data( tag, *ents_to_set, tag_storage_data );MB_CHK_ERR( rval );

    return moab::MB_SUCCESS;  // no error
}
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 2067 of file iMOAB.cpp.

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

Referenced by main().

{
#ifdef MOAB_HAVE_MPI
    appData& data = context.appDatas[*pid];
    Range ent_exchange;
    std::vector< Tag > tags;

    for( int i = 0; i < *num_tag; i++ )
    {
        if( tag_indices[i] < 0 || tag_indices[i] >= (int)data.tagList.size() )
        {
            return moab::MB_FAILURE;
        }  // error in tag index

        tags.push_back( data.tagList[tag_indices[i]] );
    }

    if( *ent_type == 0 )
    {
        ent_exchange = data.all_verts;
    }
    else if( *ent_type == 1 )
    {
        ent_exchange = data.primary_elems;
    }
    else
    {
        return moab::MB_FAILURE;
    }  // unexpected type

    ParallelComm* pco = context.pcomms[*pid];

    ErrorCode rval = pco->exchange_tags( tags, tags, ent_exchange );MB_CHK_ERR( rval );

#else
    /* do nothing if serial */
    // just silence the warning
    // do not call sync tags in serial!
    int k = *pid + *num_tag + *tag_indices + *ent_type;
    k++;
#endif

    return moab::MB_SUCCESS;
}

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: Not 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 769 of file iMOAB.cpp.

References appData::all_verts, GlobalContext::appDatas, moab::Range::clear(), 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, 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::point_cloud, appData::primary_elems, PSTATUS_GHOST, PSTATUS_NOT, moab::Range::size(), moab::subtract(), and moab::Interface::UNION.

Referenced by iMOAB_LoadMesh().

{
    // this will include ghost elements info
    appData& data        = context.appDatas[*pid];
    EntityHandle fileSet = data.file_set;

    // first clear all data ranges; this can be called after ghosting
    data.all_verts.clear();
    data.primary_elems.clear();
    data.local_verts.clear();
    data.ghost_vertices.clear();
    data.owned_elems.clear();
    data.ghost_elems.clear();
    data.mat_sets.clear();
    data.neu_sets.clear();
    data.diri_sets.clear();

    // Let us get all the vertex entities
    ErrorCode rval = context.MBI->get_entities_by_type( fileSet, MBVERTEX, data.all_verts, true );MB_CHK_ERR( rval );  // recursive

    // Let us check first entities of dimension = 3
    data.dimension = 3;
    rval           = context.MBI->get_entities_by_dimension( fileSet, data.dimension, data.primary_elems, true );MB_CHK_ERR( rval );  // recursive

    if( data.primary_elems.empty() )
    {
        // Now 3-D elements. Let us check entities of dimension = 2
        data.dimension = 2;
        rval           = context.MBI->get_entities_by_dimension( fileSet, data.dimension, data.primary_elems, true );MB_CHK_ERR( rval );  // recursive

        if( data.primary_elems.empty() )
        {
            // Now 3-D/2-D elements. Let us check entities of dimension = 1
            data.dimension = 1;
            rval = context.MBI->get_entities_by_dimension( fileSet, data.dimension, data.primary_elems, true );MB_CHK_ERR( rval );  // recursive

            if( data.primary_elems.empty() )
            {
                // no elements of dimension 1 or 2 or 3; it could happen for point clouds
                data.dimension = 0;
            }
        }
    }

    // check if the current mesh is just a point cloud
    data.point_cloud = ( ( data.primary_elems.size() == 0 && data.all_verts.size() > 0 ) || data.dimension == 0 );

#ifdef MOAB_HAVE_MPI

    if( context.MPI_initialized )
    {
        ParallelComm* pco = context.pcomms[*pid];

        // filter ghost vertices, from local
        rval = pco->filter_pstatus( data.all_verts, PSTATUS_GHOST, PSTATUS_NOT, -1, &data.local_verts );MB_CHK_ERR( rval );

        // Store handles for all ghosted entities
        data.ghost_vertices = subtract( data.all_verts, data.local_verts );

        // filter ghost elements, from local
        rval = pco->filter_pstatus( data.primary_elems, PSTATUS_GHOST, PSTATUS_NOT, -1, &data.owned_elems );MB_CHK_ERR( rval );

        data.ghost_elems = subtract( data.primary_elems, data.owned_elems );
    }
    else
    {
        data.local_verts = data.all_verts;
        data.owned_elems = data.primary_elems;
    }

#else

    data.local_verts = data.all_verts;
    data.owned_elems = data.primary_elems;

#endif

    // Get the references for some standard internal tags such as material blocks, BCs, etc
    rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.material_tag ), 0, 1,
                                                      data.mat_sets, Interface::UNION );MB_CHK_ERR( rval );

    rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.neumann_tag ), 0, 1,
                                                      data.neu_sets, Interface::UNION );MB_CHK_ERR( rval );

    rval = context.MBI->get_entities_by_type_and_tag( fileSet, MBENTITYSET, &( context.dirichlet_tag ), 0, 1,
                                                      data.diri_sets, Interface::UNION );MB_CHK_ERR( rval );

    return moab::MB_SUCCESS;
}

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 751 of file iMOAB.cpp.

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

Referenced by main().

{
    IMOAB_CHECKPOINTER( prefix, 2 );
    IMOAB_ASSERT( strlen( prefix ), "Invalid prefix string length." );

    std::ostringstream file_name;
    int rank = 0, size = 1;
#ifdef MOAB_HAVE_MPI
    rank = context.pcomms[*pid]->rank();
    size = context.pcomms[*pid]->size();
#endif
    file_name << prefix << "_" << size << "_" << rank << ".h5m";
    // Now let us actually write the file to disk with appropriate options
    ErrorCode rval = context.MBI->write_file( file_name.str().c_str(), 0, 0, &context.appDatas[*pid].file_set, 1 );MB_CHK_ERR( rval );

    return moab::MB_SUCCESS;
}
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 690 of file iMOAB.cpp.

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().

Referenced by commgraphtest(), main(), migrate(), and migrate_smart().

{
    IMOAB_CHECKPOINTER( filename, 2 );
    IMOAB_ASSERT( strlen( filename ), "Invalid filename length." );

    appData& data        = context.appDatas[*pid];
    EntityHandle fileSet = data.file_set;

    std::ostringstream newopts;
#ifdef MOAB_HAVE_MPI
    std::string write_opts( ( write_options ? write_options : "" ) );
    std::string pcid( "PARALLEL_COMM=" );
    std::size_t found = write_opts.find( pcid );

    if( found != std::string::npos )
    {
        std::cerr << " cannot specify PARALLEL_COMM option, it is implicit \n";
        return moab::MB_FAILURE;
    }

    // if write in parallel, add pc option, to be sure about which ParallelComm instance is used
    std::string pw( "PARALLEL=WRITE_PART" );
    found = write_opts.find( pw );

    if( found != std::string::npos )
    {
        newopts << "PARALLEL_COMM=" << *pid << ";";
    }

#endif

    if( write_options ) newopts << write_options;

    std::vector< Tag > copyTagList = data.tagList;
    // append Global ID and Parallel Partition
    std::string gid_name_tag( "GLOBAL_ID" );

    // export global id tag, we need it always
    if( data.tagMap.find( gid_name_tag ) == data.tagMap.end() )
    {
        Tag gid = context.MBI->globalId_tag();
        copyTagList.push_back( gid );
    }
    // also Parallel_Partition PARALLEL_PARTITION
    std::string pp_name_tag( "PARALLEL_PARTITION" );

    // write parallel part tag too, if it exists
    if( data.tagMap.find( pp_name_tag ) == data.tagMap.end() )
    {
        Tag ptag;
        context.MBI->tag_get_handle( pp_name_tag.c_str(), ptag );
        if( ptag ) copyTagList.push_back( ptag );
    }

    // Now let us actually write the file to disk with appropriate options
    ErrorCode rval = context.MBI->write_file( filename, 0, newopts.str().c_str(), &fileSet, 1, &copyTagList[0],
                                              (int)copyTagList.size() );MB_CHK_ERR( rval );

    return moab::MB_SUCCESS;
}
static void split_tag_names ( std::string  input_names,
std::string &  separator,
std::vector< std::string > &  list_tag_names 
) [static]

Definition at line 438 of file iMOAB.cpp.

Referenced by iMOAB_DefineTagStorage(), iMOAB_GetDoubleTagStorage(), iMOAB_SetDoubleTagStorage(), and iMOAB_SetDoubleTagStorageWithGid().

{
    size_t pos = 0;
    std::string token;
    while( ( pos = input_names.find( separator ) ) != std::string::npos )
    {
        token = input_names.substr( 0, pos );
        if( !token.empty() ) list_tag_names.push_back( token );
        // std::cout << token << std::endl;
        input_names.erase( 0, pos + separator.length() );
    }
    if( !input_names.empty() )
    {
        // if leftover something, or if not ended with delimiter
        list_tag_names.push_back( input_names );
    }
    return;
}

Variable Documentation

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines