Skip to content

Commit

Permalink
Fix bug in writing VTK cell data - believe will close #776
Browse files Browse the repository at this point in the history
  • Loading branch information
rupertnash committed Oct 6, 2022
1 parent 7f4b2f5 commit 3a3c7a3
Show file tree
Hide file tree
Showing 3 changed files with 51 additions and 55 deletions.
50 changes: 25 additions & 25 deletions Code/redblood/MeshIO.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,15 +18,16 @@
#include "redblood/CellEnergy.h"
#include "redblood/VTKError.h"
#include "util/fileutils.h"
#include "util/Iterator.h"

namespace hemelb {
namespace redblood {

auto MeshIO::readFile(std::string const &filename, bool fixFacetOrientation) const -> MeshPtr {
return read(Mode::file, filename, fixFacetOrientation);
return read(Storage::file, filename, fixFacetOrientation);
}
auto MeshIO::readString(std::string const &data, bool fixFacetOrientation) const -> MeshPtr {
return read(Mode::string, data, fixFacetOrientation);
return read(Storage::string, data, fixFacetOrientation);
}

void MeshIO::writeFile(std::string const &filename,
Expand All @@ -37,15 +38,15 @@ namespace hemelb {
void MeshIO::writeFile(std::string const &filename,
MeshData::Vertices const& vertices, MeshData::Facets const& facets,
util::UnitConverter const& c, PointScalarData const& data) const {
write(Mode::file, filename, vertices, facets, c, data);
write(Storage::file, filename, vertices, facets, c, data);
}

std::string MeshIO::writeString(MeshData const& mesh, util::UnitConverter const& c, PointScalarData const& data) const {
return writeString(mesh.vertices, mesh.facets, c, data);
}

std::string MeshIO::writeString(MeshData::Vertices const& vertices, MeshData::Facets const& facets, util::UnitConverter const& c, PointScalarData const& data) const {
return write(Mode::string, std::string{}, vertices, facets, c, data);
return write(Storage::string, std::string{}, vertices, facets, c, data);
}

void MeshIO::writeFile(std::string const& filename, CellBase const& cell, util::UnitConverter const& c) const {
Expand Down Expand Up @@ -191,17 +192,17 @@ namespace hemelb {
return result;
}

auto KruegerMeshIO::read(Mode mode, std::string const &filename_or_data, bool fixFacetOrientation) const -> MeshPtr {
auto KruegerMeshIO::read(Storage mode, std::string const &filename_or_data, bool fixFacetOrientation) const -> MeshPtr {
switch (mode) {
case Mode::file:
case Storage::file:
if (!util::file_exists(filename_or_data.c_str()))
throw Exception() << "Red-blood-cell mesh file '" << filename_or_data << "' does not exist";

log::Logger::Log<log::Debug, log::Singleton>("Reading red blood cell from %s",
filename_or_data.c_str());
// Open file if it exists
return read_krueger_mesh(std::ifstream{filename_or_data}, fixFacetOrientation);
case Mode::string:
case Storage::string:
return read_krueger_mesh(std::istringstream{filename_or_data}, fixFacetOrientation);
}
}
Expand Down Expand Up @@ -233,9 +234,9 @@ namespace hemelb {
stream << "$EndElement\n";
}

std::string KruegerMeshIO::write(Mode m, std::string const& filename,
MeshData::Vertices const& vertices, MeshData::Facets const& facets,
util::UnitConverter const& c, PointScalarData const& data) const {
std::string KruegerMeshIO::write(Storage m, std::string const& filename,
MeshData::Vertices const& vertices, MeshData::Facets const& facets,
util::UnitConverter const& c, PointScalarData const& data) const {
if (data.size()) {
std::string msg{"Krueger mesh IO does not support point data. Omitting fields:"};
for (auto&& pair: data) {
Expand All @@ -245,7 +246,7 @@ namespace hemelb {
}

switch (m) {
case Mode::file:
case Storage::file:
{
log::Logger::Log<log::Debug, log::Singleton>("Writing red blood cell to %s",
filename.c_str());
Expand All @@ -254,7 +255,7 @@ namespace hemelb {
write_krueger_mesh(file, vertices, facets, c);
return {};
}
case Mode::string:
case Storage::string:
{
std::ostringstream ss;
write_krueger_mesh(ss, vertices, facets, c);
Expand All @@ -268,17 +269,17 @@ namespace hemelb {
// VTK format
//

auto VTKMeshIO::readUnoriented(Mode mode, std::string const &filename_or_data) const -> std::tuple<MeshPtr, PolyDataPtr>
auto VTKMeshIO::readUnoriented(Storage mode, std::string const &filename_or_data) const -> std::tuple<MeshPtr, PolyDataPtr>
{
// Read in VTK polydata object
auto reader = ErrThrower<vtkXMLPolyDataReader>::New();
switch (mode) {
case Mode::file:
case Storage::file:
log::Logger::Log<log::Debug, log::Singleton>("Reading red blood cell from VTK polydata file");
reader->ReadFromInputStringOff();
reader->SetFileName(filename_or_data.c_str());
break;
case Mode::string:
case Storage::string:
log::Logger::Log<log::Debug, log::Singleton>("Reading red blood cell from VTK polydata string");
reader->ReadFromInputStringOn();
reader->SetInputString(filename_or_data);
Expand Down Expand Up @@ -336,7 +337,7 @@ namespace hemelb {
return std::make_tuple(result, polydata);
}

auto VTKMeshIO::read(Mode mode, std::string const &filename_or_data, bool fixFacetOrientation) const -> MeshPtr {
auto VTKMeshIO::read(Storage mode, std::string const &filename_or_data, bool fixFacetOrientation) const -> MeshPtr {
MeshIO::MeshPtr meshData;
VTKMeshIO::PolyDataPtr polyData;
std::tie(meshData, polyData) = readUnoriented(mode, filename_or_data);
Expand All @@ -349,19 +350,18 @@ namespace hemelb {
return meshData;
}

std::string VTKMeshIO::write(Mode m, std::string const &filename,
MeshData::Vertices const& vertices, MeshData::Facets const& facets,
util::UnitConverter const& c, PointScalarData const& pt_scalar_fields) const {
std::string VTKMeshIO::write(Storage m, std::string const &filename,
MeshData::Vertices const& vertices, MeshData::Facets const& facets,
util::UnitConverter const& c, PointScalarData const& pt_scalar_fields) const {
// Build the vtkPolyData
auto pd = vtkSmartPointer<vtkPolyData>::New();

// First, the points/vertices
auto points = vtkSmartPointer<vtkPoints>::New();
points->SetNumberOfPoints(vertices.size());
vtkIdType i = 0;
for (auto&& v_lat: vertices) {
auto const v = c.ConvertPositionToPhysicalUnits(v_lat);
points->SetPoint(i, v.x, v.y, v.z);
for (auto&& [i, v_lat]: util::enumerate(vertices)) {
auto const v = c.ConvertPositionToPhysicalUnits(v_lat);
points->SetPoint(i, v.x, v.y, v.z);
}
pd->SetPoints(points);

Expand Down Expand Up @@ -395,15 +395,15 @@ namespace hemelb {
// Based on the type of write, configure the writer, run it and
// return any output string.
switch (m) {
case Mode::file:
case Storage::file:
log::Logger::Log<log::Debug, log::Singleton>("Writing red blood cell to %s",
filename.c_str());
writer->WriteToOutputStringOff();
writer->SetFileName(filename.c_str());
writer->Write();
return {};

case Mode::string:
case Storage::string:
writer->WriteToOutputStringOn();
writer->Write();
return writer->GetOutputString();
Expand Down
48 changes: 24 additions & 24 deletions Code/redblood/MeshIO.h
Original file line number Diff line number Diff line change
Expand Up @@ -37,7 +37,7 @@ namespace hemelb {
using PointScalarData = std::vector<ScalarField>;

// Where are we writing to?
enum class Mode {file, string};
enum class Storage {file, string};

virtual ~MeshIO() = default;

Expand Down Expand Up @@ -79,48 +79,48 @@ namespace hemelb {

private:
// Derived class implement these to actually do the serialisation.
virtual MeshPtr read(Mode,
std::string const& filename_or_data,
bool fixFacetOrientation) const = 0;
virtual std::string write(Mode m, std::string const& filename,
MeshData::Vertices const &, MeshData::Facets const &,
util::UnitConverter const& c, PointScalarData const& data) const = 0;
virtual MeshPtr read(Storage,
std::string const& filename_or_data,
bool fixFacetOrientation) const = 0;
virtual std::string write(Storage m, std::string const& filename,
MeshData::Vertices const &, MeshData::Facets const &,
util::UnitConverter const& c, PointScalarData const& data) const = 0;
};

// Format is from T. Krueger's thesis
class KruegerMeshIO : public MeshIO {
public:
~KruegerMeshIO() = default;

MeshPtr read(Mode,
std::string const &filename_or_data,
bool fixFacetOrientation) const override;
MeshPtr read(Storage,
std::string const &filename_or_data,
bool fixFacetOrientation) const override;

std::string write(Mode m,
std::string const &filename,
MeshData::Vertices const &, MeshData::Facets const &,
util::UnitConverter const& c,
PointScalarData const& data) const override;
std::string write(Storage m,
std::string const &filename,
MeshData::Vertices const &, MeshData::Facets const &,
util::UnitConverter const& c,
PointScalarData const& data) const override;
};

// VTK XML PolyData format
class VTKMeshIO : public MeshIO {
public:
using PolyDataPtr = vtkSmartPointer<vtkPolyData>;
// Helper member function, exposed to allow testing of the fix orientations
std::tuple<MeshPtr, PolyDataPtr> readUnoriented(Mode, std::string const &) const;
std::tuple<MeshPtr, PolyDataPtr> readUnoriented(Storage, std::string const &) const;

~VTKMeshIO() = default;

MeshPtr read(Mode,
std::string const &filename_or_data,
bool fixFacetOrientation) const override;
MeshPtr read(Storage,
std::string const &filename_or_data,
bool fixFacetOrientation) const override;

std::string write(Mode m,
std::string const &filename,
MeshData::Vertices const &, MeshData::Facets const &,
util::UnitConverter const& c,
PointScalarData const& data) const override;
std::string write(Storage m,
std::string const &filename,
MeshData::Vertices const &, MeshData::Facets const &,
util::UnitConverter const& c,
PointScalarData const& data) const override;

};

Expand Down
8 changes: 2 additions & 6 deletions Code/tests/redblood/RedBloodMeshVTKDataIOTests.cc
Original file line number Diff line number Diff line change
Expand Up @@ -56,9 +56,7 @@ namespace hemelb
// This file has 684 out of 720 faces oriented inwards
std::string filename = resources::Resource("rbc_ico_720.vtp").Path();

std::shared_ptr<MeshData> meshData;
vtkSmartPointer<vtkPolyData> polyData;
std::tie(meshData, polyData) = io.readUnoriented(redblood::VTKMeshIO::Mode::file, filename);
auto [meshData, polyData] = io.readUnoriented(redblood::VTKMeshIO::Storage::file, filename);

auto numSwaps = orientFacets(*meshData, *polyData);

Expand All @@ -70,9 +68,7 @@ namespace hemelb
// This file has all 720 faces oriented inwards
std::string filename = resources::Resource("992Particles_rank3_26_t992.vtp").Path();

std::shared_ptr<MeshData> meshData;
vtkSmartPointer<vtkPolyData> polyData;
std::tie(meshData, polyData) = io.readUnoriented(redblood::VTKMeshIO::Mode::file, filename);
auto [meshData, polyData] = io.readUnoriented(redblood::VTKMeshIO::Storage::file, filename);

auto numSwaps = orientFacets(*meshData, *polyData);

Expand Down

0 comments on commit 3a3c7a3

Please sign in to comment.