MeshKit
1.0
|
00001 00008 // Test file for parallel mesh generation 00009 // It reads geometry in parallel with CGM and meshed by camal in parallel. 00010 #include <iostream> 00011 #include <sstream> 00012 00013 #include "meshkit/MKCore.hpp" 00014 #include "meshkit/ParallelMesher.hpp" 00015 #include "meshkit/SizingFunction.hpp" 00016 #include "meshkit/ModelEnt.hpp" 00017 00018 using namespace MeshKit; 00019 00020 00021 00022 #ifdef HAVE_ACIS 00023 #define DEFAULT_TEST_FILE "bricks_3.sat" 00024 #elif defined (HAVE_OCC) 00025 #define DEFAULT_TEST_FILE "bricks_3.occ" 00026 #endif 00027 00028 int load_and_mesh(const char *geom_filename, 00029 const char *mesh_filename, 00030 const char *options, const double interval_size, 00031 const int n_interval, const int rank); 00032 00033 int main( int argc, char *argv[] ) 00034 { 00035 int rank, size; 00036 MPI_Init(&argc, &argv); 00037 MPI_Comm_rank(MPI_COMM_WORLD, &rank); 00038 MPI_Comm_size(MPI_COMM_WORLD, &size); 00039 std::cout << "Hello" << std::endl; 00040 00041 // Check command line arg 00042 std::string geom_filename; 00043 const char *mesh_filename = NULL; 00044 int send_method = 1; // read and delete 00045 int part_method = 0; // round_robin 00046 double mesh_size = 0.1; 00047 int mesh_interval = -1; 00048 bool force_intervals = false; 00049 std::string options; 00050 00051 if (argc < 2) { 00052 if (rank == 0) { 00053 std::cout << "Usage: " << argv[0] << " <input_geom_filename> <output_mesh_filename> [<send_methond>] [<partition_method>] [-s <uniform_size>] [-i <uniform_int>] [-f] " 00054 << std::endl 00055 << " <send_method> = 0:read, 1:read_delete, 2:broadcast" 00056 << " , 3:bcast_delete, 4:scatter," 00057 << " , 5:read_parallel" << std::endl; 00058 std::cout << " <partition_method> = 0:round-robin, 1:static" << std::endl; 00059 std::cout << " -s <uniform_size> = mesh with this size" << std::endl; 00060 std::cout << " -i <uniform_int> = mesh curves with this # intervals" << std::endl; 00061 std::cout << " -f = force these size/interval settings even if geometry has interval settings" << std::endl; 00062 std::cout << "No file specified. Defaulting to: " << DEFAULT_TEST_FILE << std::endl; 00063 std::cout << "No send method specified. Defaulting to: read_delete" << std::endl; 00064 } 00065 std::string file_name = (std::string) MESH_DIR + "/" + DEFAULT_TEST_FILE; 00066 geom_filename += (std::string) MESH_DIR; 00067 geom_filename += "/"; 00068 geom_filename += DEFAULT_TEST_FILE; 00069 } 00070 else { 00071 geom_filename = argv[1]; 00072 if (argc > 2) mesh_filename = argv[2]; 00073 if (argc > 3) send_method = atoi(argv[3]); 00074 if (argc > 4) part_method = atoi(argv[4]); 00075 if (argc > 5) { 00076 int argno = 5; 00077 while (argno < argc) { 00078 if (!strcmp(argv[argno], "-s")) { 00079 argno++; 00080 sscanf(argv[argno], "%lf", &mesh_size); 00081 argno++; 00082 } 00083 else if (!strcmp(argv[argno], "-i")) { 00084 argno++; 00085 sscanf(argv[argno], "%d", &mesh_interval); 00086 argno++; 00087 } 00088 else if (!strcmp(argv[argno], "-f")) { 00089 argno++; 00090 force_intervals = true; 00091 } 00092 else { 00093 std::cerr << "Unrecognized option: " << argv[argno] << std::endl; 00094 return 1; 00095 } 00096 } 00097 } 00098 #ifdef HAVE_ACIS 00099 if (send_method != 1) { 00100 std::cerr << "Other send methods than read_and_delete for ACIS geometry are not supported." << std::endl; 00101 return 0; 00102 } 00103 #endif 00104 } 00105 00106 // set geometry send method 00107 if (send_method == 0) options += "PARALLEL=READ;"; 00108 else if (send_method == 1) options += "PARALLEL=READ_DELETE;"; 00109 else if (send_method == 2) options += "PARALLEL=BCAST;"; 00110 else if (send_method == 3) options += "PARALLEL=BCAST_DELETE;"; 00111 else if (send_method == 4) options += "PARALLEL=SCATTER;"; 00112 else if (send_method == 5) options += "PARALLEL=READ_PARALLEL;"; 00113 else { 00114 std::cout << "Send method " << send_method 00115 << " is not supported. Defaulting to: read_delete" << std::endl; 00116 options = "PARALLEL=READ_DELETE;"; 00117 } 00118 00119 // set partition method 00120 if (part_method == 0) options += "PARTITION=GEOM_DIMENSION;PARTITION_VAL=3;PARTITION_DISTRIBUTE;"; 00121 else if (part_method == 1) options += "PARTITION=PAR_PARTITION_STATIC;"; 00122 else { 00123 std::cout << "Partition method " << part_method 00124 << " is not supported. Defaulting to: round_robin" << std::endl; 00125 options = "PARTITION=GEOM_DIMENSION;PARTITION_VAL=3;PARTITION_DISTRIBUTE;"; 00126 } 00127 00128 if (load_and_mesh(geom_filename.c_str(), mesh_filename, 00129 options.c_str(), mesh_size, mesh_interval, rank)) return 1; 00130 00131 return 0; 00132 } 00133 00134 int load_and_mesh(const char *geom_filename, 00135 const char *mesh_filename, 00136 const char *options, const double interval_size, 00137 const int n_interval, const int rank) 00138 { 00139 // start up MK and load the geometry 00140 double t1 = MPI_Wtime(); 00141 MKCore *mk = new MKCore; 00142 mk->load_geometry(geom_filename, options); 00143 MPI_Barrier(MPI_COMM_WORLD); 00144 double t2 = MPI_Wtime(); 00145 00146 // get the volumes 00147 MEntVector dum, vols, part_vols; 00148 mk->get_entities_by_dimension(3, vols); 00149 00150 // make a sizing function and set it on the surface 00151 SizingFunction esize(mk, n_interval, interval_size); 00152 unsigned int i_sf = esize.core_index(); 00153 for (int i = 0; i < (int) vols.size(); i++) vols[i]->sizing_function_index(i_sf); 00154 00155 // do parallel mesh 00156 mk->construct_meshop("ParallelMesher", vols); 00157 mk->setup_and_execute(); 00158 double t3 = MPI_Wtime(); 00159 00160 if (mesh_filename != NULL) { 00161 std::string out_name; 00162 std::stringstream ss; 00163 ss << mesh_filename << ".h5m"; 00164 ss >> out_name; 00165 iMesh::Error err = mk->imesh_instance()->save(0, out_name.c_str(), 00166 "moab:PARALLEL=WRITE_PART"); 00167 IBERRCHK(err, "Couldn't save mesh."); 00168 double t4 = MPI_Wtime(); 00169 std::cout << "Export_time=" 00170 << t4 - t3 << std::endl; 00171 } 00172 00173 MPI_Barrier(MPI_COMM_WORLD); 00174 double t5 = MPI_Wtime(); 00175 00176 // report the number of tets 00177 moab::Range tets; 00178 mk->moab_instance()->get_entities_by_dimension(0, 3, tets); 00179 00180 std::cout << tets.size() << " tets generated." << std::endl; 00181 00182 if (rank == 0) { 00183 std::cout << "Geometry_loading_time=" << t2 - t1 00184 << " secs, Meshing_time=" << t3 - t2 00185 << " secs, Total_time=" << t5 - t1 << std::endl; 00186 } 00187 00188 mk->clear_graph(); 00189 MPI_Finalize(); 00190 00191 return 0; 00192 }