Skip to content
Snippets Groups Projects
Commit e5f8ea19 authored by Ansgar Burchardt's avatar Ansgar Burchardt
Browse files

Merge branch '67-vtkcheck-call-python-only-once' into 'master'

VTK tests: check all VTK files with one python invocation

Closes #67

See merge request !203
parents 2b362ea4 a0bebf65
No related branches found
No related tags found
1 merge request!203VTK tests: check all VTK files with one python invocation
Pipeline #
......@@ -8,12 +8,18 @@
#include <ios>
#include <iostream>
#include <ostream>
#include <set>
#include <sstream>
#include <string>
#include <sys/wait.h>
#include <dune/common/exceptions.hh>
#include <dune/common/test/testsuite.hh>
namespace Dune {
namespace Impl {
// quote so the result can be used inside '...' in python
// quotes not included in the result
......@@ -83,9 +89,9 @@ inline bool is_suffix(const std::string &s, const std::string &suffix)
s.compare(s.size() - suffix.size(), suffix.size(), suffix) == 0;
}
inline int checkVTKFile(const std::string &name)
inline bool havePythonVTK()
{
static const bool havePythonVTK = [] {
static const bool result = [] {
// This check is invoked only once, even in a multithreading environment,
// since it is invoked in the initializer of a static variable.
if(runPython("from vtk import *") == 0)
......@@ -95,32 +101,78 @@ inline int checkVTKFile(const std::string &name)
<< "vtk can read the files we wrote." << std::endl;
return false;
} ();
if(!havePythonVTK)
{
std::cerr << "skip: " << name << std::endl;
return 77;
}
std::string reader;
if (is_suffix(name, ".vtu")) reader = "vtkXMLUnstructuredGridReader";
else if(is_suffix(name, ".pvtu")) reader = "vtkXMLPUnstructuredGridReader";
else if(is_suffix(name, ".vtp")) reader = "vtkXMLPolyDataReader";
else if(is_suffix(name, ".pvtp")) reader = "vtkXMLPPolyDataReader";
else DUNE_THROW(Dune::NotImplemented,
"Unknown vtk file extension: " << name);
std::cout << "Loading " << name << " using python vtk" << std::endl;
std::string pycode =
"from vtk import *;"
"import sys;"
"reader = "+reader+"();"
"reader.SetFileName('"+pyq(name)+"');"
"reader.Update();"
// check that the number of of cells is > 0
"sys.exit(not (reader.GetOutput().GetNumberOfCells() > 0));";
int result = runPython(pycode);
std::cout << (result == 0 ? "pass: " : "fail: ") << name << std::endl;
return result;
}
inline std::string pythonVTKReader(const std::string& filename)
{
if (is_suffix(filename, ".vtu")) return "vtkXMLUnstructuredGridReader";
else if(is_suffix(filename, ".pvtu")) return "vtkXMLPUnstructuredGridReader";
else if(is_suffix(filename, ".vtp")) return "vtkXMLPolyDataReader";
else if(is_suffix(filename, ".pvtp")) return "vtkXMLPPolyDataReader";
else DUNE_THROW(Dune::NotImplemented,
"Unknown vtk file extension: " << filename);
}
} /* namespace Impl */
class VTKChecker
{
public:
void push(const std::string& file)
{
auto res = files_.insert(file);
if (not res.second) {
testSuite_.check(false, "VTKChecker")
<< "'" << file << "' was added multiple times";
}
}
int check()
{
if (not Impl::havePythonVTK()) {
return 77;
}
else if (not files_.empty()) {
const int result = Impl::runPython(generatePythonCode());
testSuite_.check(result == 0);
}
return testSuite_.exit();
}
const TestSuite& testSuite() const
{
return testSuite_;
}
private:
std::string generatePythonCode() const
{
std::stringstream code;
code << "from vtk import *\n"
<< "import sys\n"
<< "passed = True\n";
for (const auto& file : files_) {
code << "reader = " << Impl::pythonVTKReader(file) << "()\n"
<< "reader.SetFileName('" << Impl::pyq(file) << "')\n"
<< "reader.Update()\n"
<< "if (not (reader.GetOutput().GetNumberOfCells() > 0)):\n"
<< " print('ERROR in {}'.format('" << Impl::pyq(file) << "'))\n"
<< " passed = False\n";
}
code << "sys.exit(0 if passed else 1)\n";
return code.str();
}
std::set< std::string > files_;
TestSuite testSuite_;
};
} /* namespace Dune */
#endif // DUNE_GRID_IO_FILE_TEST_CHECKVTKFILE_HH
......@@ -84,7 +84,7 @@ struct Acc
};
template< class GridView >
int doWrite( const GridView &gridView, bool coerceToSimplex)
int doWrite( Dune::VTKChecker& vtkChecker, const std::string& gridViewName, const GridView &gridView, bool coerceToSimplex)
{
enum { dim = GridView :: dimension };
......@@ -104,21 +104,21 @@ int doWrite( const GridView &gridView, bool coerceToSimplex)
int result = 0;
std::string name;
std::ostringstream prefix;
prefix << "subsamplingvtktest-" << dim << "D-"
prefix << "subsamplingvtktest-" << dim << "D-" << gridViewName << "-"
<< (coerceToSimplex ? "simplex" : "natural");
int rank = gridView.comm().rank();
name = vtk.write(prefix.str() + "-ascii");
if(rank == 0) acc(result, checkVTKFile(name));
if(rank == 0) vtkChecker.push(name);
name = vtk.write(prefix.str() + "-appendedraw", Dune::VTK::appendedraw);
if(rank == 0) acc(result, checkVTKFile(name));
if(rank == 0) vtkChecker.push(name);
return result;
}
template<int dim>
int vtkCheck(const std::array<int, dim>& elements,
int vtkCheck(Dune::VTKChecker& vtkChecker, const std::array<int, dim>& elements,
const Dune::FieldVector<double, dim>& upperRight)
{
Dune::YaspGrid<dim> g(upperRight, elements);
......@@ -132,13 +132,13 @@ int vtkCheck(const std::array<int, dim>& elements,
int result = 0;
acc(result, doWrite( g.leafGridView(), false));
acc(result, doWrite( g.levelGridView( 0 ), false));
acc(result, doWrite( g.levelGridView( g.maxLevel() ), false));
acc(result, doWrite( vtkChecker, "leafview", g.leafGridView(), false));
acc(result, doWrite( vtkChecker, "coarselevelview", g.levelGridView( 0 ), false));
acc(result, doWrite( vtkChecker, "finelevelview", g.levelGridView( g.maxLevel() ), false));
acc(result, doWrite( g.leafGridView(), true));
acc(result, doWrite( g.levelGridView( 0 ), true));
acc(result, doWrite( g.levelGridView( g.maxLevel() ), true));
acc(result, doWrite( vtkChecker, "leafview", g.leafGridView(), true));
acc(result, doWrite( vtkChecker, "coarselevelview", g.levelGridView( 0 ), true));
acc(result, doWrite( vtkChecker, "finelevelview", g.levelGridView( g.maxLevel() ), true));
return result;
}
......@@ -156,9 +156,13 @@ int main(int argc, char **argv)
int result = 0; // pass by default
using Dune::make_array;
acc(result, vtkCheck<1>(make_array(5), {1.0}));
acc(result, vtkCheck<2>(make_array(5,5), {1.0, 2.0}));
acc(result, vtkCheck<3>(make_array(5,5,5), {1.0, 2.0, 3.0}));
Dune::VTKChecker vtkChecker;
acc(result, vtkCheck<1>(vtkChecker, make_array(5), {1.0}));
acc(result, vtkCheck<2>(vtkChecker, make_array(5,5), {1.0, 2.0}));
acc(result, vtkCheck<3>(vtkChecker, make_array(5,5,5), {1.0, 2.0, 3.0}));
acc(result, vtkChecker.check());
mpiHelper.getCollectiveCommunication().allreduce<Acc>(&result, 1);
......
......@@ -94,7 +94,7 @@ struct Acc
};
template< class GridView >
int doWrite( const GridView &gridView, Dune :: VTK :: DataMode dm )
int doWrite( Dune::VTKChecker& vtkChecker, const std::string& gridViewName, const GridView &gridView, Dune :: VTK :: DataMode dm )
{
enum { dim = GridView :: dimension };
......@@ -113,27 +113,27 @@ int doWrite( const GridView &gridView, Dune :: VTK :: DataMode dm )
int result = 0;
std::string name;
std::ostringstream prefix;
prefix << "vtktest-" << dim << "D-" << VTKDataMode(dm);
prefix << "vtktest-" << dim << "D-" << gridViewName << "-" << VTKDataMode(dm);
int rank = gridView.comm().rank();
name = vtk.write(prefix.str() + "-ascii");
if(rank == 0) acc(result, checkVTKFile(name));
if(rank == 0) vtkChecker.push(name);
name = vtk.write(prefix.str() + "-base64", Dune::VTK::base64);
if(rank == 0) acc(result, checkVTKFile(name));
if(rank == 0) vtkChecker.push(name);
name = vtk.write(prefix.str() + "-appendedraw", Dune::VTK::appendedraw);
if(rank == 0) acc(result, checkVTKFile(name));
if(rank == 0) vtkChecker.push(name);
name = vtk.write(prefix.str() + "-appendedbase64",
Dune::VTK::appendedbase64);
if(rank == 0) acc(result, checkVTKFile(name));
if(rank == 0) vtkChecker.push(name);
return result;
}
template<int dim>
int vtkCheck(const std::array<int, dim>& elements,
int vtkCheck(Dune::VTKChecker& vtkChecker, const std::array<int, dim>& elements,
const Dune::FieldVector<double, dim>& upperRight)
{
Dune::YaspGrid<dim> g(upperRight, elements);
......@@ -147,15 +147,15 @@ int vtkCheck(const std::array<int, dim>& elements,
int result = 0;
acc(result, doWrite( g.leafGridView(), Dune::VTK::conforming ));
acc(result, doWrite( g.leafGridView(), Dune::VTK::nonconforming ));
acc(result, doWrite( g.levelGridView( 0 ),
acc(result, doWrite( vtkChecker, "leafview", g.leafGridView(), Dune::VTK::conforming ));
acc(result, doWrite( vtkChecker, "leafview", g.leafGridView(), Dune::VTK::nonconforming ));
acc(result, doWrite( vtkChecker, "coarselevelview", g.levelGridView( 0 ),
Dune::VTK::conforming ));
acc(result, doWrite( g.levelGridView( 0 ),
acc(result, doWrite( vtkChecker, "coarselevelview", g.levelGridView( 0 ),
Dune::VTK::nonconforming ));
acc(result, doWrite( g.levelGridView( g.maxLevel() ),
acc(result, doWrite( vtkChecker, "finelevelview", g.levelGridView( g.maxLevel() ),
Dune::VTK::conforming ));
acc(result, doWrite( g.levelGridView( g.maxLevel() ),
acc(result, doWrite( vtkChecker, "finelevelview", g.levelGridView( g.maxLevel() ),
Dune::VTK::nonconforming ));
return result;
......@@ -174,9 +174,13 @@ int main(int argc, char **argv)
int result = 0; // pass by default
using Dune::make_array;
acc(result, vtkCheck<1>(make_array(5), {1.0}));
acc(result, vtkCheck<2>(make_array(5,5), {1.0, 2.0}));
acc(result, vtkCheck<3>(make_array(5,5,5), {1.0, 2.0, 3.0}));
Dune::VTKChecker vtkChecker;
acc(result, vtkCheck<1>(vtkChecker, make_array(5), {1.0}));
acc(result, vtkCheck<2>(vtkChecker, make_array(5,5), {1.0, 2.0}));
acc(result, vtkCheck<3>(vtkChecker, make_array(5,5,5), {1.0, 2.0, 3.0}));
acc(result, vtkChecker.check());
mpiHelper.getCollectiveCommunication().allreduce<Acc>(&result, 1);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment