MOAB: Mesh Oriented datABase
(version 5.4.1)
|
00001 #include "TestUtil.hpp" 00002 #include "moab/Core.hpp" 00003 #include "moab/ParallelComm.hpp" 00004 #include "moab/ProgOptions.hpp" 00005 #include "MBParallelConventions.h" 00006 #include "moab/ReadUtilIface.hpp" 00007 00008 using namespace moab; 00009 00010 std::string example = TestDir + "unittest/io/homme3x3458.t.3.nc"; 00011 00012 void test_read_parallel_ucd_trivial(); 00013 void test_read_parallel_ucd_trivial_spectral(); 00014 void test_read_parallel( int num_verts, bool test_nb_nodes ); 00015 00016 void test_multiple_loads_of_same_file(); 00017 00018 std::string partition_method; 00019 const int levels = 3; 00020 00021 int main( int argc, char* argv[] ) 00022 { 00023 MPI_Init( &argc, &argv ); 00024 int result = 0; 00025 00026 result += RUN_TEST( test_read_parallel_ucd_trivial ); 00027 result += RUN_TEST( test_read_parallel_ucd_trivial_spectral ); 00028 result += RUN_TEST( test_multiple_loads_of_same_file ); 00029 00030 MPI_Finalize(); 00031 return result; 00032 } 00033 00034 void test_read_parallel_ucd_trivial() 00035 { 00036 partition_method = std::string( ";PARTITION_METHOD=TRIVIAL;PARALLEL_RESOLVE_SHARED_ENTS" ); 00037 test_read_parallel( 3458, true ); 00038 } 00039 00040 void test_read_parallel_ucd_trivial_spectral() 00041 { 00042 partition_method = std::string( ";PARTITION_METHOD=TRIVIAL;SPECTRAL_MESH;PARALLEL_RESOLVE_SHARED_ENTS" ); 00043 test_read_parallel( 3458, false ); 00044 } 00045 00046 void test_read_parallel( int num_verts, bool test_nb_nodes ) 00047 { 00048 Core moab; 00049 Interface& mb = moab; 00050 EntityHandle file_set; 00051 ErrorCode rval; 00052 rval = mb.create_meshset( MESHSET_SET, file_set );CHECK_ERR( rval ); 00053 00054 std::string opt = std::string( "PARALLEL=READ_PART" ) + partition_method; 00055 // Create gather set in processor 0 00056 opt += std::string( ";GATHER_SET=0" ); 00057 rval = mb.load_file( example.c_str(), &file_set, opt.c_str() );CHECK_ERR( rval ); 00058 00059 ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 ); 00060 int procs = pcomm->proc_config().proc_size(); 00061 int rank = pcomm->proc_config().proc_rank(); 00062 00063 rval = pcomm->check_all_shared_handles();CHECK_ERR( rval ); 00064 00065 // Get the total # owned verts 00066 Range verts; 00067 rval = mb.get_entities_by_type( 0, MBVERTEX, verts );CHECK_ERR( rval ); 00068 00069 int my_num = verts.size(); 00070 if( test_nb_nodes && 2 == procs ) 00071 { 00072 if( 0 == rank ) 00073 CHECK_EQUAL( 5283, my_num ); // Gather set vertices included 00074 else if( 1 == rank ) 00075 CHECK_EQUAL( 1825, my_num ); // Not owned vertices included 00076 } 00077 00078 rval = pcomm->filter_pstatus( verts, PSTATUS_NOT_OWNED, PSTATUS_NOT );CHECK_ERR( rval ); 00079 00080 my_num = verts.size(); 00081 if( test_nb_nodes && 2 == procs ) 00082 { 00083 if( 0 == rank ) 00084 CHECK_EQUAL( 5283, my_num ); // Gather set vertices included 00085 else if( 1 == rank ) 00086 CHECK_EQUAL( 1633, my_num ); // Not owned vertices excluded 00087 } 00088 00089 if( 0 == rank ) 00090 { 00091 // Get gather set 00092 EntityHandle gather_set; 00093 ReadUtilIface* readUtilIface; 00094 rval = mb.query_interface( readUtilIface );CHECK_ERR( rval ); 00095 rval = readUtilIface->get_gather_set( gather_set );CHECK_ERR( rval ); 00096 00097 // Get gather set entities 00098 Range gather_ents; 00099 rval = mb.get_entities_by_handle( gather_set, gather_ents );CHECK_ERR( rval ); 00100 00101 // Remove gather set vertices in processor 0 00102 verts = subtract( verts, gather_ents ); 00103 } 00104 00105 my_num = verts.size(); 00106 if( test_nb_nodes && 2 == procs ) 00107 { 00108 if( 0 == rank ) 00109 CHECK_EQUAL( 1825, my_num ); // Gather set vertices excluded 00110 else if( 1 == rank ) 00111 CHECK_EQUAL( 1633, my_num ); // Not owned vertices excluded 00112 } 00113 00114 std::cout << "proc: " << rank << " verts:" << my_num << "\n"; 00115 00116 int total_verts; 00117 MPI_Reduce( &my_num, &total_verts, 1, MPI_INT, MPI_SUM, 0, pcomm->proc_config().proc_comm() ); 00118 00119 if( 0 == rank ) 00120 { 00121 std::cout << "total vertices: " << total_verts << "\n"; 00122 if( test_nb_nodes ) CHECK_EQUAL( total_verts, num_verts ); 00123 } 00124 00125 #ifdef MOAB_HAVE_HDF5_PARALLEL 00126 std::string write_options( "PARALLEL=WRITE_PART;" ); 00127 mb.write_file( "test.h5m", NULL, write_options.c_str() ); 00128 #endif 00129 } 00130 00131 void test_multiple_loads_of_same_file() 00132 { 00133 Core moab; 00134 Interface& mb = moab; 00135 EntityHandle file_set; 00136 ErrorCode rval; 00137 rval = mb.create_meshset( MESHSET_SET, file_set );CHECK_ERR( rval ); 00138 00139 // Read first only header information, no mesh, no variable 00140 std::string opts( "PARALLEL=READ_PART;PARTITION;NOMESH;VARIABLE=;PARTITION_METHOD=TRIVIAL" ); 00141 rval = mb.load_file( example.c_str(), &file_set, opts.c_str() );CHECK_ERR( rval ); 00142 00143 opts = "PARALLEL=READ_PART;PARTITION;PARALLEL_RESOLVE_SHARED_ENTS;PARTITION_METHOD=TRIVIAL;" 00144 "VARIABLE="; 00145 // Create gather set in processor 1 00146 opts += std::string( ";GATHER_SET=1" ); 00147 rval = mb.load_file( example.c_str(), &file_set, opts.c_str() );CHECK_ERR( rval ); 00148 00149 opts = "PARALLEL=READ_PART;PARTITION;PARTITION_METHOD=TRIVIAL;NOMESH;VARIABLE=T;TIMESTEP=0"; 00150 rval = mb.load_file( example.c_str(), &file_set, opts.c_str() );CHECK_ERR( rval ); 00151 00152 // Check values of tag T0 (first level) at some strategically chosen places below 00153 ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 ); 00154 int procs = pcomm->proc_config().proc_size(); 00155 int rank = pcomm->proc_config().proc_rank(); 00156 00157 // Make check runs this test in two processors 00158 if( 2 == procs ) 00159 { 00160 Range verts; 00161 rval = mb.get_entities_by_type( 0, MBVERTEX, verts );CHECK_ERR( rval ); 00162 00163 int my_num = verts.size(); 00164 if( 0 == rank ) 00165 CHECK_EQUAL( 1825, my_num ); 00166 else if( 1 == rank ) 00167 { 00168 CHECK_EQUAL( 5283, 00169 my_num ); // Gather set vertices included; Not owned vertices included 00170 00171 // Get gather set 00172 EntityHandle gather_set; 00173 ReadUtilIface* readUtilIface; 00174 rval = mb.query_interface( readUtilIface );CHECK_ERR( rval ); 00175 rval = readUtilIface->get_gather_set( gather_set );CHECK_ERR( rval ); 00176 00177 // Get gather set entities 00178 Range gather_ents; 00179 rval = mb.get_entities_by_handle( gather_set, gather_ents );CHECK_ERR( rval ); 00180 00181 // Remove gather set vertices in processor 1 00182 verts = subtract( verts, gather_ents ); 00183 my_num = verts.size(); 00184 CHECK_EQUAL( 1825, 00185 my_num ); // Gather set vertices excluded; Not owned vertices included 00186 } 00187 00188 Tag Ttag0; 00189 rval = mb.tag_get_handle( "T0", levels, MB_TYPE_DOUBLE, Ttag0, MB_TAG_DENSE );CHECK_ERR( rval ); 00190 00191 int count; 00192 void* Tbuf; 00193 rval = mb.tag_iterate( Ttag0, verts.begin(), verts.end(), count, Tbuf );CHECK_ERR( rval ); 00194 CHECK_EQUAL( count, my_num ); 00195 00196 const double eps = 0.0001; 00197 double* data = (double*)Tbuf; 00198 00199 if( 0 == rank ) 00200 { 00201 CHECK_REAL_EQUAL( 233.1136, data[0 * levels], eps ); // First vert 00202 CHECK_REAL_EQUAL( 237.1977, data[912 * levels], eps ); // Median vert 00203 CHECK_REAL_EQUAL( 234.9711, data[1824 * levels], eps ); // Last vert 00204 } 00205 else if( 1 == rank ) 00206 { 00207 CHECK_REAL_EQUAL( 233.1136, data[0 * levels], eps ); // First vert 00208 CHECK_REAL_EQUAL( 231.0446, data[912 * levels], eps ); // Median vert 00209 CHECK_REAL_EQUAL( 234.0416, data[1824 * levels], eps ); // Last vert 00210 } 00211 } 00212 }