test/parallel-export.cpp

00001 #include "DAPTA.hpp"
00002 
00003 #include <set>
00004 #include <fstream>
00005 #include <iostream>
00006 #include <iomanip>
00007 #undef SEEK_SET
00008 #undef SEEK_CUR
00009 #undef SEEK_END
00010 #include "mpi.h"
00011 using namespace DAPTA;
00012 
00013 #define LENGTH_EXPORT_FILE 0
00014 #define STRING_EXPORT_FILE 1
00015 
00016 typedef Mesh<Tri3MeshTraits> TriMesh;
00017 TriMesh mesh;
00018 RCBDecomposer<TriMesh::VertexType> rcb;
00019 std::vector<TriMesh::VertexHandle> points;
00020 
00021 int main(int argc, char **argv) {
00022    int          myrank, nprocs;
00023    MPI_Status   status;
00024    
00025    // Initialize MPI
00026    MPI_Init(&argc, &argv);
00027 
00028    // Get my rank id and number of MPI processes
00029    MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
00030    MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
00031    
00032    if(myrank == 0) {
00033       // Load in the mesh
00034       mesh.ClearMesh();
00035       std::ifstream meshfile;
00036       int num;
00037       
00038       meshfile.open(argv[1]);
00039       if(!meshfile) {
00040          std::cout << "ERROR READING FILE!" << std::endl;
00041          exit(1);
00042       }
00043       
00044       mesh.ReadFile(meshfile);
00045       meshfile.close();
00046       points = mesh.GetCentres();
00047       
00048       rcb.loadVertices(points);
00049       
00050       for(int i=1; i<nprocs; i++) {
00051          //std::cout << "client(" << myrank << "): GETTING DOMAIN FOR DOMAIN: " << i << "/" << nprocs << std::endl;
00052          
00053          std::vector<TriMesh::VertexHandle> myNodes = rcb.decompose(nprocs, i);
00054          std::set<TriMesh::VertexHandle> part_vert;
00055          std::set<TriMesh::ElementHandle> part_elem;
00056          std::string export_file;
00057          
00058          for(int j=0; j<myNodes.size(); j++) {
00059             //std::cout << "   client(" << myrank << "): DOING ELEMENT NUMBER: " << myNodes.at(j)->globalID << std::endl;
00060             
00061             TriMesh::ElementHandle this_elem = mesh.FindElementFromIndex(myNodes.at(j)->globalID);
00062             part_elem.insert(this_elem);
00063             for(int k=0; k<TriMesh::ElementType::num_vertices; k++) {
00064                part_vert.insert(this_elem->vertices.at(k));
00065             }
00066          }
00067          
00068          //std::cout << "client(" << myrank << "): GETTING EXPORT FILE FOR DOMAIN: " << i << "/" << nprocs << std::endl;         
00069          export_file = mesh.ExportDomainFile(part_vert, part_elem);
00070 
00071          const char* char_send = export_file.c_str();
00072          //strcpy(char_send, export_file.c_str());
00073          int length_export_file = strlen(char_send);
00074          
00075          //std::cout << "client(" << myrank << "): SENDING LENGTH OF EXPORT FILE FOR DOMAIN: " << i << "/" << nprocs << std::endl;
00076          //std::cout << "   client(" << myrank << "): LENGTH IS: " << length_export_file << std::endl;
00077          MPI_Send(&length_export_file, 1, MPI_INT, i, LENGTH_EXPORT_FILE, MPI_COMM_WORLD);
00078 
00079          //std::cout << "client(" << myrank << "): SENDING EXPORT FILE FOR DOMAIN: " << i << "/" << nprocs << std::endl;
00080          MPI_Send((void *)char_send, length_export_file, MPI_CHAR, i, STRING_EXPORT_FILE, MPI_COMM_WORLD);
00081          
00082          //std::cout << "client(" << myrank << "): DONE SENDING TO DOMAIN: " << i << "/" << nprocs << std::endl;
00083       }
00084       
00085       //std::cout << "client(" << myrank << "): NOW LOADING IN THE MESH FILE..." << std::endl;
00086       std::vector<TriMesh::VertexHandle> myNodes = rcb.decompose(nprocs, 0);
00087       myNodes = rcb.decompose(nprocs, 0);
00088       std::set<TriMesh::VertexHandle> part_vert;
00089       std::set<TriMesh::ElementHandle> part_elem;
00090       std::string export_file;
00091       
00092       for(int j=0; j<myNodes.size(); j++) {
00093          TriMesh::ElementHandle this_elem = mesh.FindElementFromIndex(myNodes.at(j)->globalID);
00094          part_elem.insert(this_elem);
00095          for(int k=0; k<TriMesh::ElementType::num_vertices; k++) {
00096             part_vert.insert(this_elem->vertices.at(k));
00097          }
00098       }
00099       
00100       //std::cout << "client(" << myrank << "): READING MESH..." << std::endl;
00101       export_file = mesh.ExportDomainFile(part_vert, part_elem);
00102       std::stringstream export_file_stream;
00103       export_file_stream.str(export_file);
00104       
00105       mesh.ClearMesh();
00106       mesh.ReadFile(export_file_stream);
00107    }
00108    else {
00109       //std::cout << "client(" << myrank << "): RECEIVING DATA..." << std::endl;
00110       int  length_export_file;
00111       char export_file_char[] = "";
00112       //std::cout << "client(" << myrank << "): RECEIVING LENGTH" << std::endl;
00113       MPI_Recv(&length_export_file,                1,  MPI_INT, 0, LENGTH_EXPORT_FILE, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
00114       //std::cout << "   client(" << myrank << "): RECEIVED LENGTH: " << length_export_file << std::endl;
00115       //std::cout << "client(" << myrank << "): RECEIVING EXPORT FILE" << std::endl;
00116       MPI_Recv(export_file_char, length_export_file, MPI_CHAR, 0, STRING_EXPORT_FILE, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
00117       //std::cout << "client(" << myrank << "): RECEIVED EXPORT FILE" << std::endl;
00118       
00119       //std::cout << "client(" << myrank << "): LOADING THE MESH..." << std::endl;
00120       std::string       export_file(export_file_char);
00121       std::stringstream export_file_stream;
00122       export_file_stream.str(export_file);
00123       
00124       //std::cout << "client(" << myrank << "): CLEARING MESH..." << std::endl;
00125       mesh.ClearMesh();
00126       //std::cout << "client(" << myrank << "): READING MESH..." << std::endl;
00127       mesh.ReadFile(export_file_stream);
00128    }
00129    
00130    std::cout << "client(" << myrank << "): ================================================" << std::endl;
00131    for(TriMesh::VertexMapIterator v_iter=mesh.VertexMapBegin(); v_iter!=mesh.VertexMapEnd(); ++v_iter) {
00132       std::cout << std::setprecision(4);
00133       std::cout << "client(" << myrank << "): " << (*v_iter)->globalID << "\t";
00134       for(int i=0; i<TriMesh::VertexType::dimension; i++) {
00135          std::cout << (*v_iter)->coords.at(i) << "\t";
00136       }
00137       std::cout << (*v_iter) << std::endl;
00138    }
00139    std::cout << "client(" << myrank << "): " << std::endl;
00140    for(TriMesh::ElementMapIterator e_iter=mesh.ElementMapBegin(); e_iter!=mesh.ElementMapEnd(); ++e_iter) {
00141       std::cout << "client(" << myrank << "): " << (*e_iter)->globalID << "\t";
00142       for(int i=0; i<TriMesh::ElementType::num_vertices; i++) {
00143          std::cout << (*e_iter)->vertices.at(i)->globalID << "\t";
00144       }
00145       std::cout << (*e_iter) << std::endl;
00146    }
00147    std::cout << "client(" << myrank << "): ================================================" << std::endl;
00148    
00149    //std::cout << "client(" << myrank << "): WAITING FOR ALL PROCESSES..." << std::endl;
00150    MPI_Barrier(MPI_COMM_WORLD);
00151    
00152    // Shut Down MPI
00153    MPI_Finalize();
00154    
00155    return 0;
00156 }

Generated on Tue May 29 17:13:49 2007 for DAPTA by  doxygen 1.5.1