easy-iso

Iso-surface extraction from volume data
git clone https://git.0xfab.ch/easy-iso.git
Log | Files | Refs | Submodules | README | LICENSE

commit 97a49114f5c8e123c6d9363e34379c21c1d83124
parent ed7dc0c0333eed97a8ff17b31336e299e77409d2
Author: Fabian Wermelinger <fabianw@mavt.ethz.ch>
Date:   Thu,  1 Dec 2016 12:22:53 +0100

added parallel MPI I/O

Diffstat:
Msrc/IsoExtractor.h | 6+++---
Msrc/IsoExtractorMPI.h | 117+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++----------------
2 files changed, 97 insertions(+), 26 deletions(-)

diff --git a/src/IsoExtractor.h b/src/IsoExtractor.h @@ -18,8 +18,8 @@ struct Triangle { - Real vertices[3][3]; -}; + float vertices[3][3]; +} __attribute__((packed)); static const Real gVertexOffset[8][3] = { @@ -532,7 +532,7 @@ void IsoExtractor<Tkernel,Tinterpolator>::_extractIso(const Real isoval, std::chrono::time_point<std::chrono::system_clock> start, end; start = std::chrono::system_clock::now(); if (verbose) - std::cout << "Extracting isovalue " << isoval << "..." << std::endl; + std::cout << "Extracting isovalue " << isoval << std::endl; const size_t nCubes = static_cast<size_t>(Nx) * Ny * Nz; size_t totalSize = 0; diff --git a/src/IsoExtractorMPI.h b/src/IsoExtractorMPI.h @@ -10,6 +10,12 @@ #include "IsoExtractor.h" #include "InterpolatorMPI.h" +struct TriangleConnect +{ + unsigned char n; + int v[3]; +} __attribute__((packed)); + template <typename Tkernel, template<typename> class Tinterpolator=InterpolatorMPI> class IsoExtractorMPI : public IsoExtractor<Tkernel,Tinterpolator> { @@ -109,7 +115,8 @@ private: // write content // _writePLY_binary_MPI(comm, tList); - _writePLY_binary_MPI_lazy(tList, basename+".ply"); + // _writePLY_binary_MPI_lazy(tList, basename+".ply"); + _writePLY_binary_MPI(comm, tList, basename+".ply"); } void _writePLY_binary_MPI_lazy(const std::vector<Triangle>& tList, const std::string fname) @@ -164,28 +171,92 @@ private: } } - // void _writePLY_binary_MPI(const MPI_Comm comm, const std::vector<Triangle>& tList) - // { - // for (auto t : tList) - // { - // float tc[3][3]; - // for (int j = 0; j < 3; ++j) - // { - // tc[j][0] = t.vertices[j][0]; - // tc[j][1] = t.vertices[j][1]; - // tc[j][2] = t.vertices[j][2]; - // } - // ply.write((char*)&tc[0][0], 9*sizeof(float)); - // } - - // const char nVert = 3; - // for (int i = 0; i < tList.size()*3; i += 3) - // { - // int tv[3] = {i, i+1, i+2}; - // ply.write(&nVert, sizeof(char)); - // ply.write((char*)&tv[0], 3*sizeof(int)); - // } - // } + void _writePLY_binary_MPI(const MPI_Comm comm, const std::vector<Triangle>& tList, const std::string fname) + { + int rank, size; + MPI_Status status; + MPI_File file_id; + MPI_Comm_rank(comm, &rank); + MPI_Comm_size(comm, &size); + + int rc = MPI_File_open(MPI_COMM_SELF, fname.c_str(), MPI_MODE_APPEND | MPI_MODE_WRONLY, MPI_INFO_NULL, &file_id); + if (rc) + { + if (0 == rank) + std::cerr << "ERROR: IsoExtractorMPI: Can not open file." << std::endl; + abort(); + } + + std::vector<unsigned long long> allsize(size); + unsigned long long mytri = tList.size(); + MPI_Allgather(&mytri, 1, MPI_UNSIGNED_LONG_LONG, allsize.data(), 1, MPI_UNSIGNED_LONG_LONG, comm); + + // if (0 == rank) + // { + // for (const auto v: allsize) + // std::cout << v << std::endl; + // } + + long long offset = 0; + MPI_File_get_position(file_id, &offset); + + assert(sizeof(Triangle) == 9*sizeof(float)); + for (int i = 0; i < rank; ++i) + offset = offset + allsize[i]*sizeof(Triangle); + + Triangle* tarr = nullptr; + if (allsize[rank] > 0) + { + tarr = new Triangle[allsize[rank]]; + // for (size_t i=0; i<allsize[rank]; ++i) + // memcpy(&tarr[0].vertices[0][0] + i*sizeof(Triangle), (char*)&tList[i].vertices[0][0], sizeof(Triangle)); + for (size_t i=0; i<allsize[rank]; ++i) + { + for (int j = 0; j < 3; ++j) + { + tarr[i].vertices[j][0] = tList[i].vertices[j][0]; + tarr[i].vertices[j][1] = tList[i].vertices[j][1]; + tarr[i].vertices[j][2] = tList[i].vertices[j][2]; + } + } + + MPI_File_write_at(file_id, offset, (char *)&tarr[0].vertices[0][0], allsize[rank]*sizeof(Triangle), MPI_CHAR, &status); + } + offset = offset + allsize[rank]*sizeof(Triangle); + MPI_Barrier(comm); + MPI_Bcast(&offset, 1, MPI_LONG_LONG, size-1, comm); + + if (tarr != nullptr) + delete [] tarr; + + for (int i = 0; i < rank; ++i) + offset = offset + allsize[i]*sizeof(TriangleConnect); + + TriangleConnect* tcarr = nullptr; + if (allsize[rank] > 0) + { + tcarr = new TriangleConnect[allsize[rank]]; + int base = 0; + for (int i = 0; i < rank; ++i) + base += 3*allsize[i]; + + for (size_t i=0; i<allsize[rank]; ++i) + { + tcarr[i].n = 3; + tcarr[i].v[0] = base + 3*i + 0; + tcarr[i].v[1] = base + 3*i + 1; + tcarr[i].v[2] = base + 3*i + 2; + } + + MPI_File_write_at(file_id, offset, (char *)&tcarr[0].n, allsize[rank]*sizeof(TriangleConnect), MPI_CHAR, &status); + } + MPI_Barrier(comm); + + if (tcarr != nullptr) + delete [] tcarr; + + MPI_File_close(&file_id); + } }; #endif /* ISOEXTRACTORMPI_H_XABTB0WN */