MOAB: Mesh Oriented datABase  (version 5.4.1)
migrate_nontrivial.cpp File Reference
#include "moab/ParallelComm.hpp"
#include "moab/Core.hpp"
#include "moab_mpi.h"
#include "moab/iMOAB.h"
#include "TestUtil.hpp"
#include "moab/ProgOptions.hpp"
+ Include dependency graph for migrate_nontrivial.cpp:

Go to the source code of this file.

Defines

#define RUN_TEST_ARG2(A, B)   run_test( &( A ), #A, B )
#define CHECKRC(rc, message)

Functions

int is_any_proc_error (int is_my_error)
int run_test (ErrorCode(*func)(const char *), const char *func_name, const char *file_name)
ErrorCode migrate_graph (const char *filename)
ErrorCode migrate_geom (const char *filename)
ErrorCode migrate_trivial (const char *filename)
ErrorCode migrate_smart (const char *filename, const char *outfile, int partMethod)
int main (int argc, char *argv[])

Variables

int rank
int size
int ierr
int compid1
int compid2
int nghlay
std::vector< int > groupTasks
int startG1
int startG2
int endG1
int endG2
MPI_Comm jcomm
MPI_Group jgroup

Define Documentation

#define CHECKRC (   rc,
  message 
)
Value:
if( 0 != ( rc ) )                     \
    {                                     \
        printf( "Error: %s\n", message ); \
        return MB_FAILURE;                \
    }

Definition at line 30 of file migrate_nontrivial.cpp.

Referenced by migrate_smart().

#define RUN_TEST_ARG2 (   A,
 
)    run_test( &( A ), #A, B )

Definition at line 24 of file migrate_nontrivial.cpp.

Referenced by main().


Function Documentation

int is_any_proc_error ( int  is_my_error)

Definition at line 37 of file migrate_nontrivial.cpp.

References MPI_COMM_WORLD.

{
    int result = 0;
    int err    = MPI_Allreduce( &is_my_error, &result, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD );
    return err || result;
}
int main ( int  argc,
char *  argv[] 
)

Definition at line 207 of file migrate_nontrivial.cpp.

References ProgOptions::addOpt(), endG1, endG2, filename, jcomm, jgroup, migrate_geom(), migrate_graph(), migrate_trivial(), MPI_COMM_WORLD, ProgOptions::parseCommandLine(), rank, RUN_TEST_ARG2, size, startG1, and startG2.

{
    MPI_Init( &argc, &argv );
    MPI_Comm_rank( MPI_COMM_WORLD, &rank );
    MPI_Comm_size( MPI_COMM_WORLD, &size );

    MPI_Comm_dup( MPI_COMM_WORLD, &jcomm );
    MPI_Comm_group( jcomm, &jgroup );

    ProgOptions opts;
    int typeTest = 3;
    // std::string inputfile, outfile("out.h5m"), netcdfFile, variable_name, sefile_name;
    std::string filename;
    filename = TestDir + "unittest/field1.h5m";
    startG1  = 0;
    startG2  = 0;
    endG1    = 0;
    endG2    = 1;

    opts.addOpt< std::string >( "file,f", "source file", &filename );

    opts.addOpt< int >( "startSender,a", "start task for source layout", &startG1 );
    opts.addOpt< int >( "endSender,b", "end task for source layout", &endG1 );
    opts.addOpt< int >( "startRecv,c", "start task for receiver layout", &startG2 );
    opts.addOpt< int >( "endRecv,d", "end task for receiver layout", &endG2 );

    opts.addOpt< int >( "typeTest,t", "test types (0 - trivial, 1 graph, 2 geom, 3 both  graph and geometry",
                        &typeTest );

    opts.parseCommandLine( argc, argv );

    if( rank == 0 )
    {
        std::cout << " input file : " << filename << "\n";
        std::cout << " sender   on tasks: " << startG1 << ":" << endG1 << "\n";
        std::cout << " receiver on tasks: " << startG2 << ":" << endG2 << "\n";
        std::cout << " type migrate: " << typeTest << " (0 - trivial, 1 graph , 2 geom, 3 both graph and geom  ) \n";
    }

    int num_errors = 0;

    if( 0 == typeTest ) num_errors += RUN_TEST_ARG2( migrate_trivial, filename.c_str() );
    if( 3 == typeTest || 1 == typeTest ) num_errors += RUN_TEST_ARG2( migrate_graph, filename.c_str() );
    if( 3 == typeTest || 2 == typeTest ) num_errors += RUN_TEST_ARG2( migrate_geom, filename.c_str() );

    if( rank == 0 )
    {
        if( !num_errors )
            std::cout << "All tests passed" << std::endl;
        else
            std::cout << num_errors << " TESTS FAILED!" << std::endl;
    }

    MPI_Group_free( &jgroup );
    MPI_Comm_free( &jcomm );
    MPI_Finalize();
    return num_errors;
}
ErrorCode migrate_geom ( const char *  filename)

Definition at line 197 of file migrate_nontrivial.cpp.

References migrate_smart().

Referenced by main().

{
    return migrate_smart( filename, "migrate_geom.h5m", 2 );
}
ErrorCode migrate_graph ( const char *  filename)

Definition at line 192 of file migrate_nontrivial.cpp.

References migrate_smart().

Referenced by main().

{
    return migrate_smart( filename, "migrate_graph.h5m", 1 );
}
ErrorCode migrate_smart ( const char *  filename,
const char *  outfile,
int  partMethod 
)

Definition at line 76 of file migrate_nontrivial.cpp.

References CHECKRC, compid1, compid2, context, endG1, endG2, groupTasks, ierr, iMOAB_AppID, iMOAB_DeregisterApplication(), iMOAB_Finalize(), iMOAB_Initialize(), iMOAB_LoadMesh(), iMOAB_RegisterApplication(), iMOAB_WriteMesh(), jcomm, jgroup, MB_SUCCESS, nghlay, readopts(), startG1, and startG2.

Referenced by migrate_geom(), migrate_graph(), and migrate_trivial().

{
    // first create MPI groups

    std::string filen( filename );
    MPI_Group group1, group2;
    groupTasks.resize( endG1 - startG1 + 1 );
    for( int i = startG1; i <= endG1; i++ )
        groupTasks[i - startG1] = i;

    ierr = MPI_Group_incl( jgroup, endG1 - startG1 + 1, &groupTasks[0], &group1 );
    CHECKRC( ierr, "can't create group1" )

    groupTasks.resize( endG2 - startG2 + 1 );
    for( int i = startG2; i <= endG2; i++ )
        groupTasks[i - startG2] = i;

    ierr = MPI_Group_incl( jgroup, endG2 - startG2 + 1, &groupTasks[0], &group2 );
    CHECKRC( ierr, "can't create group2" )

    // create 2 communicators, one for each group
    int tagcomm1 = 1, tagcomm2 = 2;
    int context_id = -1;  // plain migrate, default context
    MPI_Comm comm1, comm2;
    ierr = MPI_Comm_create_group( jcomm, group1, tagcomm1, &comm1 );
    CHECKRC( ierr, "can't create comm1" )

    ierr = MPI_Comm_create_group( jcomm, group2, tagcomm2, &comm2 );
    CHECKRC( ierr, "can't create comm2" )

    ierr = iMOAB_Initialize( 0, 0 );  // not really needed anything from argc, argv, yet; maybe we should
    CHECKRC( ierr, "can't initialize iMOAB" )

    // give some dummy values to component ids, just to differentiate between them
    // the par comm graph is unique between components
    compid1 = 4;
    compid2 = 7;

    int appID1;
    iMOAB_AppID pid1 = &appID1;
    int appID2;
    iMOAB_AppID pid2 = &appID2;

    if( comm1 != MPI_COMM_NULL )
    {
        ierr = iMOAB_RegisterApplication( "APP1", &comm1, &compid1, pid1 );
        CHECKRC( ierr, "can't register app1 " )
    }
    if( comm2 != MPI_COMM_NULL )
    {
        ierr = iMOAB_RegisterApplication( "APP2", &comm2, &compid2, pid2 );
        CHECKRC( ierr, "can't register app2 " )
    }

    if( comm1 != MPI_COMM_NULL )
    {

        std::string readopts( "PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARALLEL_RESOLVE_SHARED_ENTS" );

        nghlay = 0;

        ierr = iMOAB_LoadMesh( pid1, filen.c_str(), readopts.c_str(), &nghlay );
        CHECKRC( ierr, "can't load mesh " )
        ierr = iMOAB_SendMesh( pid1, &jcomm, &group2, &compid2, &partMethod );  // send to component 2
        CHECKRC( ierr, "cannot send elements" )
#ifdef GRAPH_INFO
        int is_sender = 1;
        int context   = compid2;
        iMOAB_DumpCommGraph( pid1, &context, &is_sender, "MigrateS" );
#endif
    }

    if( comm2 != MPI_COMM_NULL )
    {
        ierr = iMOAB_ReceiveMesh( pid2, &jcomm, &group1, &compid1 );  // receive from component 1
        CHECKRC( ierr, "cannot receive elements" )
        std::string wopts;
        wopts = "PARALLEL=WRITE_PART;";
        ierr  = iMOAB_WriteMesh( pid2, outfile, wopts.c_str() );
        CHECKRC( ierr, "cannot write received mesh" )
#ifdef GRAPH_INFO
        int is_sender = 0;
        int context   = compid1;
        iMOAB_DumpCommGraph( pid2, &context, &is_sender, "MigrateR" );
#endif
    }

    MPI_Barrier( jcomm );

    // we can now free the sender buffers
    context_id = compid2;  // even for default migrate, be more explicit
    if( comm1 != MPI_COMM_NULL ) ierr = iMOAB_FreeSenderBuffers( pid1, &context_id );

    if( comm2 != MPI_COMM_NULL )
    {
        ierr = iMOAB_DeregisterApplication( pid2 );
        CHECKRC( ierr, "cannot deregister app 2 receiver" )
    }

    if( comm1 != MPI_COMM_NULL )
    {
        ierr = iMOAB_DeregisterApplication( pid1 );
        CHECKRC( ierr, "cannot deregister app 1 sender" )
    }

    ierr = iMOAB_Finalize();
    CHECKRC( ierr, "did not finalize iMOAB" )

    if( MPI_COMM_NULL != comm1 ) MPI_Comm_free( &comm1 );
    if( MPI_COMM_NULL != comm2 ) MPI_Comm_free( &comm2 );

    MPI_Group_free( &group1 );
    MPI_Group_free( &group2 );
    return MB_SUCCESS;
}
ErrorCode migrate_trivial ( const char *  filename)

Definition at line 202 of file migrate_nontrivial.cpp.

References migrate_smart().

Referenced by main().

{
    return migrate_smart( filename, "migrate_trivial.h5m", 0 );
}
int run_test ( ErrorCode(*)(const char *)  func,
const char *  func_name,
const char *  file_name 
)

Definition at line 44 of file migrate_nontrivial.cpp.

References ErrorCode, is_any_proc_error(), MB_SUCCESS, MPI_COMM_WORLD, and rank.

{
    ErrorCode result = ( *func )( file_name );
    int is_err       = is_any_proc_error( ( MB_SUCCESS != result ) );
    int rank;
    MPI_Comm_rank( MPI_COMM_WORLD, &rank );
    if( rank == 0 )
    {
        if( is_err )
            std::cout << func_name << " : FAILED!!" << std::endl;
        else
            std::cout << func_name << " : success" << std::endl;
    }

    return is_err;
}

Variable Documentation

int compid1

Definition at line 68 of file migrate_nontrivial.cpp.

Referenced by migrate(), and migrate_smart().

int compid2

Definition at line 68 of file migrate_nontrivial.cpp.

Referenced by migrate(), and migrate_smart().

int endG1

Definition at line 71 of file migrate_nontrivial.cpp.

int endG2

Definition at line 71 of file migrate_nontrivial.cpp.

std::vector< int > groupTasks

Definition at line 70 of file migrate_nontrivial.cpp.

int ierr

Definition at line 66 of file migrate_nontrivial.cpp.

Definition at line 73 of file migrate_nontrivial.cpp.

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

MPI_Group jgroup

Definition at line 74 of file migrate_nontrivial.cpp.

int nghlay

Definition at line 69 of file migrate_nontrivial.cpp.

int rank
Examples:
CrystalRouterExample.cpp, ErrorHandlingSimulation.cpp, GenLargeMesh.cpp, HelloParMOAB.cpp, ReadWriteTest.cpp, ReduceExchangeTags.cpp, and StructuredMeshSimple.cpp.

Definition at line 66 of file migrate_nontrivial.cpp.

Referenced by moab::MeshOutputFunctor::assign_global_ids(), moab::MeshGeneration::BrickInstance(), PartMap::build_map(), check_consistent_ids(), closedsurface_uref_hirec_convergence_study(), moab::NCWriteGCRM::collect_mesh_info(), moab::NCWriteHOMME::collect_mesh_info(), moab::NCWriteMPAS::collect_mesh_info(), moab::WriteHDF5Parallel::communicate_shared_set_data(), PartMap::count_from_rank(), create_coarse_mesh(), moab::WriteHDF5Parallel::create_dataset(), moab::NCHelperScrip::create_mesh(), moab::NCHelperHOMME::create_mesh(), moab::NCHelperGCRM::create_mesh(), moab::NCHelperMPAS::create_mesh(), moab::WriteHDF5Parallel::create_meshset_tables(), create_shared_grid_3d(), moab::WriteHDF5Parallel::create_tag_tables(), moab::Coupler::do_normalization(), do_rank_subst(), moab::HiReconstruction::eval_vander_bivar_cmf(), moab::HiReconstruction::eval_vander_univar_cmf(), moab::ReadHDF5::find_sets_containing(), gather_one_cell_var(), generate_mesh(), moab::TempestOnlineMap::GenerateRemappingWeights(), GenerateTestMatrixAndVectors(), get_attrib_array_length_handle(), moab::Coupler::get_matching_entities(), get_max_id(), moab::SharedSetData::get_owner_handle(), get_tag(), handle_error_code(), iMOAB_LoadMesh(), iMOAB_WriteLocalMesh(), moab::NCHelperDomain::init_mesh_vals(), moab::NCHelperFV::init_mesh_vals(), moab::NCHelperEuler::init_mesh_vals(), initialize_tree(), interface_verts(), moab::TempestOnlineMap::IsConservative(), moab::TempestOnlineMap::LinearRemapFVtoFV_Tempest_MOAB(), moab::ReadParallel::load_file(), load_meshset_hirec(), main(), moab::MBTraceBackErrorHandler(), mhdf_getTagInfo(), mhdf_open_table(), MPI_swap(), multiple_loads_of_same_file(), moab::my_Gatherv(), parallel_create_mesh(), moab::ReadNC::parse_options(), moab::WriteNC::parse_options(), PartMap::part_from_coords(), ZoltanPartitioner::partition_owned_cells(), print_partitioned_entities(), moab::ProcConfig::ProcConfig(), moab::ReadHDF5::read_all_set_meta(), read_buffered_map(), read_mesh_parallel(), read_one_cell_var(), moab::TempestOnlineMap::ReadParallelMap(), ZoltanPartitioner::repartition(), report_iface_ents(), report_nsets(), resolve_and_exchange(), run_test(), runner_run_tests(), set_local_domain_bounds(), moab::ErrorOutput::set_rank(), moab::DebugOutput::set_rank(), moab::ReadHDF5::set_up_read(), moab::TempestOnlineMap::SetDOFmapAssociation(), moab::SharedSetData::SharedSetData(), test_assign_global_ids(), test_correct_ghost(), test_entity_copies(), test_entity_copy_parts(), test_entity_owner(), test_entity_status(), test_eul_check_across_files(), test_eul_check_append(), test_eul_check_T(), test_eul_check_timestep(), test_exchange_ents(), test_fv_check_T(), test_gather_onevar(), test_gcrm_check_vars(), test_get_neighbors(), test_get_part_boundary(), test_get_parts(), test_ghost_tag_exchange(), test_ghosted_entity_shared_data(), test_homme_check_T(), test_interface_owners_common(), test_intx_in_parallel_elem_based(), test_mpas_check_vars(), test_part_boundary_iter(), test_part_id_handle(), test_part_rank(), test_push_tag_data_common(), test_read_and_ghost_after(), test_read_bc_sets(), test_read_elements_common(), test_read_eul_onevar(), test_read_fv_onevar(), test_read_global_tags(), test_read_non_adjs_side(), test_read_parallel(), test_read_sets_common(), test_read_tags(), test_read_time(), test_read_with_ghost(), test_read_with_ghost_no_augment(), test_read_with_thin_ghost_layer(), test_reduce_tag_explicit_dest(), test_reduce_tag_failures(), test_reduce_tags(), test_shared_sets(), test_string_rank_subst(), test_tempest_map_bcast(), test_trivial_partition(), test_var_length_parallel(), test_write_dense_tags(), test_write_different_element_types(), test_write_different_tags(), test_write_elements(), test_write_shared_sets(), test_write_unbalanced(), TestMeshRefiner(), tprint(), moab::NCWriteHelper::write_set_variables(), moab::TempestOnlineMap::WriteHDF5MapFile(), and moab::TempestOnlineMap::WriteSCRIPMapFile().

int size

Definition at line 66 of file migrate_nontrivial.cpp.

int startG1

Definition at line 71 of file migrate_nontrivial.cpp.

int startG2

Definition at line 71 of file migrate_nontrivial.cpp.

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines