MOAB: Mesh Oriented datABase
(version 5.4.1)
|
00001 #include "moab/ParallelComm.hpp" 00002 #include "moab/Core.hpp" 00003 #include "moab_mpi.h" 00004 #include "TestUtil.hpp" 00005 #include "MBTagConventions.hpp" 00006 #include <iostream> 00007 #include <sstream> 00008 00009 std::string filename = TestDir + "unittest/io/p8ex1.h5m"; 00010 00011 using namespace moab; 00012 void report_sets( moab::Core* mb, int rank, int nproc ) 00013 { 00014 // check neumann and material sets, and see if their number of quads / hexes in them 00015 const char* const shared_set_tag_names[] = { MATERIAL_SET_TAG_NAME, DIRICHLET_SET_TAG_NAME, NEUMANN_SET_TAG_NAME, 00016 PARALLEL_PARTITION_TAG_NAME }; 00017 00018 int num_tags = sizeof( shared_set_tag_names ) / sizeof( shared_set_tag_names[0] ); 00019 00020 for( int p = 0; p < nproc; p++ ) 00021 { 00022 if( rank == p ) 00023 { 00024 std::cout << " Task no:" << rank << "\n"; 00025 for( int i = 0; i < num_tags; i++ ) 00026 { 00027 Tag tag; 00028 ErrorCode rval = mb->tag_get_handle( shared_set_tag_names[i], 1, MB_TYPE_INTEGER, tag, MB_TAG_ANY );CHECK_ERR( rval ); 00029 Range sets; 00030 rval = mb->get_entities_by_type_and_tag( 0, MBENTITYSET, &tag, 0, 1, sets, Interface::UNION );CHECK_ERR( rval ); 00031 00032 std::vector< int > vals( sets.size() ); 00033 rval = mb->tag_get_data( tag, sets, &vals[0] );CHECK_ERR( rval ); 00034 std::cout << " sets: " << shared_set_tag_names[i] << "\n"; 00035 int j = 0; 00036 for( Range::iterator it = sets.begin(); it != sets.end(); it++, j++ ) 00037 { 00038 Range ents; 00039 rval = mb->get_entities_by_handle( *it, ents );CHECK_ERR( rval ); 00040 std::cout << " set " << mb->id_from_handle( *it ) << " with tagval=" << vals[j] << " has " 00041 << ents.size() << " entities\n"; 00042 } 00043 } 00044 } 00045 MPI_Barrier( MPI_COMM_WORLD ); 00046 } 00047 } 00048 void test_read_with_ghost() 00049 { 00050 int nproc, rank; 00051 MPI_Comm_size( MPI_COMM_WORLD, &nproc ); 00052 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 00053 moab::Core* mb = new moab::Core(); 00054 00055 ErrorCode rval = MB_SUCCESS; 00056 00057 char read_opts[] = "PARALLEL=READ_PART;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=3.0.1.3;" 00058 "PARTITION=PARALLEL_PARTITION"; 00059 rval = mb->load_file( filename.c_str(), 0, read_opts );CHECK_ERR( rval ); 00060 00061 report_sets( mb, rank, nproc ); 00062 // write in serial the database , on each rank 00063 std::ostringstream outfile; 00064 outfile << "testReadGhost_n" << nproc << "." << rank << ".h5m"; 00065 // the mesh contains ghosts too, but they were not part of mat/neumann set 00066 // write in serial the file, to see what tags are missing / or not 00067 rval = mb->write_file( outfile.str().c_str() ); // everything on root 00068 CHECK_ERR( rval ); 00069 delete mb; 00070 } 00071 00072 void test_read_and_ghost_after() 00073 { 00074 int nproc, rank; 00075 MPI_Comm_size( MPI_COMM_WORLD, &nproc ); 00076 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 00077 moab::Core* mb = new moab::Core(); 00078 moab::ParallelComm* pc = new moab::ParallelComm( mb, MPI_COMM_WORLD ); 00079 ErrorCode rval = MB_SUCCESS; 00080 00081 // first read in parallel, then ghost, then augment 00082 char read_opts[] = "PARALLEL=READ_PART;PARALLEL_RESOLVE_SHARED_ENTS;PARTITION=PARALLEL_PARTITION"; 00083 rval = mb->load_file( filename.c_str(), 0, read_opts );CHECK_ERR( rval ); 00084 00085 int ghost_dim = 3, bridge = 0, layers = 1, addl_ents = 3; 00086 rval = pc->exchange_ghost_cells( ghost_dim, bridge, layers, addl_ents, true, true );CHECK_ERR( rval ); 00087 00088 // 00089 pc->set_debug_verbosity( 1 ); 00090 rval = pc->augment_default_sets_with_ghosts( 0 );CHECK_ERR( rval ); 00091 00092 report_sets( mb, rank, nproc ); 00093 00094 // write in serial the database , on each rank 00095 std::ostringstream outfile; 00096 outfile << "TaskMesh_n" << nproc << "." << rank << ".h5m"; 00097 // the mesh contains ghosts too, but they were not part of mat/neumann set 00098 // write in serial the file, to see what tags are missing / or not 00099 rval = mb->write_file( outfile.str().c_str() ); // everything on root 00100 CHECK_ERR( rval ); 00101 00102 delete pc; 00103 delete mb; 00104 } 00105 00106 void test_read_with_ghost_no_augment() 00107 { 00108 int nproc, rank; 00109 MPI_Comm_size( MPI_COMM_WORLD, &nproc ); 00110 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 00111 moab::Core* mb = new moab::Core(); 00112 00113 ErrorCode rval = MB_SUCCESS; 00114 00115 char read_opts[] = "PARALLEL=READ_PART;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=3.0.1.3;" 00116 "PARTITION=PARALLEL_PARTITION;SKIP_AUGMENT_WITH_GHOSTS"; 00117 rval = mb->load_file( filename.c_str(), 0, read_opts );CHECK_ERR( rval ); 00118 00119 report_sets( mb, rank, nproc ); 00120 // write in serial the database , on each rank 00121 std::ostringstream outfile; 00122 outfile << "testReadGhost_n" << nproc << "." << rank << ".h5m"; 00123 // the mesh contains ghosts too, but they were not part of mat/neumann set 00124 // write in serial the file, to see what tags are missing / or not 00125 rval = mb->write_file( outfile.str().c_str() ); // everything on root 00126 CHECK_ERR( rval ); 00127 delete mb; 00128 } 00129 int main( int argc, char* argv[] ) 00130 { 00131 MPI_Init( &argc, &argv ); 00132 00133 int result = 0; 00134 if( argc >= 2 ) filename = argv[1]; 00135 00136 result += RUN_TEST( test_read_with_ghost ); 00137 result += RUN_TEST( test_read_and_ghost_after ); 00138 result += RUN_TEST( test_read_with_ghost_no_augment ); 00139 00140 MPI_Finalize(); 00141 return 0; 00142 }