Skip to content

Commit

Permalink
Change reporting of VTK errors
Browse files Browse the repository at this point in the history
  • Loading branch information
rupertnash committed Sep 11, 2023
1 parent 856d5d3 commit 1823f89
Show file tree
Hide file tree
Showing 3 changed files with 115 additions and 139 deletions.
157 changes: 80 additions & 77 deletions Code/redblood/MeshIO.cc
Original file line number Diff line number Diff line change
Expand Up @@ -195,18 +195,18 @@ namespace hemelb::redblood {
}

auto KruegerMeshIO::read(Storage mode, std::string const &filename_or_data, bool fixFacetOrientation) const -> MeshPtr {
switch (mode) {
case Storage::file:
if (!std::filesystem::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 Storage::string:
return read_krueger_mesh(std::istringstream{filename_or_data}, fixFacetOrientation);
}
switch (mode) {
case Storage::file:
if (!std::filesystem::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 Storage::string:
return read_krueger_mesh(std::istringstream{filename_or_data}, fixFacetOrientation);
}
}

static void write_krueger_mesh(std::ostream &stream, MeshData::Vertices const& vertices, MeshData::Facets const& facets, util::UnitConverter const& converter)
Expand Down Expand Up @@ -273,20 +273,21 @@ namespace hemelb::redblood {

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 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 Storage::string:
log::Logger::Log<log::Debug, log::Singleton>("Reading red blood cell from VTK polydata string");
reader->ReadFromInputStringOn();
reader->SetInputString(filename_or_data);
break;
}
// Read in VTK polydata object
VtkErrorsThrow t;
auto reader = vtkSmartPointer<vtkXMLPolyDataReader>::New();
switch (mode) {
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 Storage::string:
log::Logger::Log<log::Debug, log::Singleton>("Reading red blood cell from VTK polydata string");
reader->ReadFromInputStringOn();
reader->SetInputString(filename_or_data);
break;
}

reader->Update();

Expand Down Expand Up @@ -353,62 +354,64 @@ namespace hemelb::redblood {
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(ssize(vertices));
for (auto&& [i, v_lat]: util::enumerate(vertices)) {
auto const v = c.ConvertPositionToPhysicalUnits(v_lat);
points->SetPoint(i, v.m_values.data());
}
pd->SetPoints(points);
VtkErrorsThrow t;

// Second, the polys/facets/triangles
auto tris = vtkSmartPointer<vtkCellArray>::New();
tris->AllocateExact(ssize(facets), 3);
for (auto&& tri: facets) {
tris->InsertNextCell({tri[0], tri[1], tri[2]});
}
pd->SetPolys(tris);

// Third, the scalar point data
for (auto&& field: pt_scalar_fields) {
auto&& name = field.first;
auto&& data = field.second;
auto da = vtkSmartPointer<vtkDoubleArray>::New();
da->SetName(name.c_str());
da->SetNumberOfComponents(1);
// This let's VTK "borrow" the data inside our vector. We
// know we're not modifying the polydata, so the const cast is
// OK.
da->SetArray(const_cast<double*>(data.data()), ssize(data), 1);
pd->GetPointData()->AddArray(da);
}
// Build the vtkPolyData
auto pd = vtkSmartPointer<vtkPolyData>::New();

// Now, the writer
auto writer = ErrThrower<vtkXMLPolyDataWriter>::New();
// The vtkPolyData we just made
writer->SetInputData(pd);
// First, the points/vertices
auto points = vtkSmartPointer<vtkPoints>::New();
points->SetNumberOfPoints(ssize(vertices));
for (auto&& [i, v_lat]: util::enumerate(vertices)) {
auto const v = c.ConvertPositionToPhysicalUnits(v_lat);
points->SetPoint(i, v.m_values.data());
}
pd->SetPoints(points);

// Based on the type of write, configure the writer, run it and
// return any output string.
switch (m) {
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 {};
// Second, the polys/facets/triangles
auto tris = vtkSmartPointer<vtkCellArray>::New();
tris->AllocateExact(ssize(facets), 3);
for (auto&& tri: facets) {
tris->InsertNextCell({tri[0], tri[1], tri[2]});
}
pd->SetPolys(tris);

// Third, the scalar point data
for (auto&& field: pt_scalar_fields) {
auto&& name = field.first;
auto&& data = field.second;
auto da = vtkSmartPointer<vtkDoubleArray>::New();
da->SetName(name.c_str());
da->SetNumberOfComponents(1);
// This let's VTK "borrow" the data inside our vector. We
// know we're not modifying the polydata, so the const cast is
// OK.
da->SetArray(const_cast<double*>(data.data()), ssize(data), 1);
pd->GetPointData()->AddArray(da);
}

case Storage::string:
writer->WriteToOutputStringOn();
writer->Write();
return writer->GetOutputString();
// Now, the writer
auto writer = vtkSmartPointer<vtkXMLPolyDataWriter>::New();
// The vtkPolyData we just made
writer->SetInputData(pd);

// Based on the type of write, configure the writer, run it and
// return any output string.
switch (m) {
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 Storage::string:
writer->WriteToOutputStringOn();
writer->Write();
return writer->GetOutputString();

}
}
}

}
30 changes: 14 additions & 16 deletions Code/redblood/VTKError.cc
Original file line number Diff line number Diff line change
Expand Up @@ -4,27 +4,25 @@
// license in the file LICENSE.

#include "redblood/VTKError.h"
#include <vtkObjectFactory.h>

namespace hemelb {
namespace redblood {
namespace hemelb::redblood {

ErrLogger* ErrLogger::New() {
return new ErrLogger;
}
vtkStandardNewMacro(ErrThrowOutputWindow);

ErrThrowOutputWindow::ErrThrowOutputWindow() = default;

void ErrLogger::Clear() {
*this = ErrLogger{};
ErrThrowOutputWindow::~ErrThrowOutputWindow() = default;

void ErrThrowOutputWindow::DisplayErrorText(char const* txt) {
throw Exception() << txt;
}

void ErrLogger::Execute(vtkObject *vtkNotUsed(caller),
unsigned long event,
void *calldata)
{
if (event == vtkCommand::ErrorEvent) {
err_occurred = true;
err_message = static_cast<const char *>(calldata);
}
VtkErrorsThrow::VtkErrorsThrow() : originalWindow(vtkOutputWindow::GetInstance()) {
vtkOutputWindow::SetInstance(ErrThrowOutputWindow::New());
}

}
VtkErrorsThrow::~VtkErrorsThrow() {
vtkOutputWindow::SetInstance(originalWindow);
}
}
67 changes: 21 additions & 46 deletions Code/redblood/VTKError.h
Original file line number Diff line number Diff line change
Expand Up @@ -8,61 +8,36 @@

#include <string>

#include <vtkCommand.h>
#include <vtkOutputWindow.h>
#include <vtkSmartPointer.h>

#include "Exception.h"

namespace hemelb {
namespace redblood {

// This class implements VTK's observer interface to be notified
// in case of error.
// Only holds the last error - can reset with Clear()
class ErrLogger : public vtkCommand {
public:
bool err_occurred = false;
std::string err_message;

static ErrLogger* New();
namespace hemelb::redblood {

void Clear();
// Subclass of the default output window that turns error
// messages (e.g. from vtkErrorMacro) into HemeLB exceptions.
class ErrThrowOutputWindow : public vtkOutputWindow
{
public:
vtkTypeMacro(ErrThrowOutputWindow, vtkOutputWindow);
static ErrThrowOutputWindow* New();

void Execute(vtkObject *vtkNotUsed(caller),
unsigned long event,
void *calldata) override;
void DisplayErrorText(const char* text) override;
protected:
ErrThrowOutputWindow();
virtual ~ErrThrowOutputWindow();
};

// Wrapper for VTK objects that installs an ErrLogger observer and
// checks for errors on use, throwing if one occurred.
//
// Mostly treat just like a vtkSmartPointer<T> but you can always
// get the members with the real pointer (val) and the ErrPointer
// (log).
template <typename T>
struct ErrThrower {
using ValType = T;
using ValPointer = vtkSmartPointer<T>;
using ErrPointer = vtkSmartPointer<ErrLogger>;

ValPointer val = nullptr;
ErrPointer log = nullptr;

// Factory
static ErrThrower New() {
ErrThrower ans{ValPointer::New(), ErrPointer::New()};
ans.val->AddObserver(vtkCommand::ErrorEvent, ans.log);
return ans;
}

// Use the wrapped object
ValPointer operator->() {
if (log->err_occurred)
throw Exception() << "VTK error occurred: " << log->err_message;
return val;
}
// RAII type that replaces the VTK global output window with an instance of the above.
class VtkErrorsThrow {
vtkSmartPointer<vtkOutputWindow> originalWindow;
public:
VtkErrorsThrow();
~VtkErrorsThrow();
VtkErrorsThrow(VtkErrorsThrow const&) = delete;
VtkErrorsThrow& operator=(VtkErrorsThrow const&) = delete;
};

}
}
#endif

0 comments on commit 1823f89

Please sign in to comment.