Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Typedefs | Functions
VtkToFld.cpp File Reference
#include <vtkPolyDataReader.h>
#include <vtkPolyData.h>
#include <vtkPointData.h>
#include <vtkPoints.h>
#include <vtkCellArray.h>
#include <vtkCellDataToPointData.h>
#include <vtkContourFilter.h>
#include <boost/unordered_set.hpp>
#include <boost/program_options.hpp>
#include <LibUtilities/BasicUtils/SessionReader.h>
#include <LibUtilities/Communication/Comm.h>
#include <SpatialDomains/MeshGraph2D.h>
#include <MultiRegions/ExpList2D.h>
#include <LocalRegions/Expansion2D.h>
Include dependency graph for VtkToFld.cpp:

Go to the source code of this file.

Classes

struct  Vertex
 Represents a vertex in the mesh. More...
struct  VertexHash
 Hash function for the Vertex struct used for defining sets. More...

Typedefs

typedef boost::shared_ptr< VertexVertexSharedPtr
typedef boost::unordered_set
< VertexSharedPtr, VertexHash
VertexSet
 Define a set of Vertices with associated hash function.

Functions

bool operator== (const VertexSharedPtr &v1, const VertexSharedPtr &v2)
 Define comparison operator for the vertex struct.
int main (int argc, char *argv[])

Typedef Documentation

typedef boost::unordered_set<VertexSharedPtr, VertexHash> VertexSet

Define a set of Vertices with associated hash function.

Definition at line 73 of file VtkToFld.cpp.

typedef boost::shared_ptr<Vertex> VertexSharedPtr

Definition at line 48 of file VtkToFld.cpp.

Function Documentation

int main ( int  argc,
char *  argv[] 
)

Main function.

Usage: VtkToFld session.xml input.vtk output.fld [options]

Definition at line 81 of file VtkToFld.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::GetCoords(), Nektar::iterator, and Nektar::LibUtilities::FieldIO::Write().

{
// Set up available options
po::options_description desc("Available options");
desc.add_options()
("help,h", "Produce this help message.")
("name,n", po::value<string>()->default_value("Intensity"),
"Name of field in VTK file to use for intensity.")
("outname,m", po::value<string>()->default_value("intensity"),
"Name of field in output FLD file.")
("precision,p", po::value<double>()->default_value(1),
"Precision of vertex matching.");
po::options_description hidden("Hidden options");
hidden.add_options()
("file", po::value<vector<string> >(), "Input filename");
po::options_description cmdline_options;
cmdline_options.add(desc).add(hidden);
po::positional_options_description p;
p.add("file", -1);
po::variables_map vm;
// Parse command-line options
try
{
po::store(po::command_line_parser(argc, argv).
options(cmdline_options).positional(p).run(), vm);
po::notify(vm);
}
catch (const std::exception& e)
{
cerr << e.what() << endl;
cerr << desc;
return 1;
}
if ( vm.count("help") || vm.count("file") == 0 ||
vm["file"].as<vector<string> >().size() != 3) {
cerr << "Usage: VtkToFld session.xml intensity.vtk output.fld [options]"
<< endl;
cerr << desc;
return 1;
}
// Extract command-line argument values
std::vector<std::string> vFiles = vm["file"].as<vector<string> >();
const string infile = vFiles[1];
const string outfile = vFiles[2];
const double factor = vm["precision"].as<double>();
const string name = vm["name"].as<string>();
const string outname = vm["outname"].as<string>();
std::vector<std::string> vFilenames;
vFilenames.push_back(vFiles[0]);
vSession = LibUtilities::SessionReader::CreateInstance(2, argv, vFilenames);
try
{
//----------------------------------------------
// Read in mesh from input file
//----------------------------------------------
//----------------------------------------------
// Define Expansion
AllocateSharedPtr(vSession,graph2D);
//----------------------------------------------
//----------------------------------------------
// Set up coordinates of mesh
int coordim = Exp->GetCoordim(0);
int nq = Exp->GetNpoints();
Array<OneD, NekDouble> xc0(nq,0.0);
Array<OneD, NekDouble> xc1(nq,0.0);
Array<OneD, NekDouble> xc2(nq,0.0);
switch(coordim)
{
case 2:
Exp->GetCoords(xc0,xc1);
break;
case 3:
Exp->GetCoords(xc0,xc1,xc2);
break;
default:
ASSERTL0(false,"Coordim not valid");
break;
}
//----------------------------------------------
vtkPolyDataReader *vtkMeshReader = vtkPolyDataReader::New();
vtkMeshReader->SetFileName(infile.c_str());
vtkMeshReader->Update();
vtkPolyData *vtkMesh = vtkMeshReader->GetOutput();
vtkCellDataToPointData* c2p = vtkCellDataToPointData::New();
#if VTK_MAJOR_VERSION <= 5
c2p->SetInput(vtkMesh);
#else
c2p->SetInputData(vtkMesh);
#endif
c2p->PassCellDataOn();
c2p->Update();
vtkPolyData *vtkDataAtPoints = c2p->GetPolyDataOutput();
vtkPoints *vtkPoints = vtkMesh->GetPoints();
ASSERTL0(vtkPoints, "ERROR: cannot get points from mesh.");
vtkCellArray *vtkPolys = vtkMesh->GetPolys();
ASSERTL0(vtkPolys, "ERROR: cannot get polygons from mesh.");
vtkPointData *vtkPData = vtkDataAtPoints->GetPointData();
ASSERTL0(vtkPolys, "ERROR: cannot get point data from file.");
VertexSet points;
double p[3];
double val;
double x, y, z;
int coeff_idx;
int i,j,n;
if (!vtkDataAtPoints->GetPointData()->HasArray(name.c_str())) {
n = vtkDataAtPoints->GetPointData()->GetNumberOfArrays();
cerr << "Input file '" << infile
<< "' does not have a field named '"
<< name << "'" << endl;
cerr << "There are " << n << " arrays in this file." << endl;
for (int i = 0; i < n; ++i)
{
cerr << " "
<< vtkDataAtPoints->GetPointData()->GetArray(i)->GetName()
<< endl;
}
return 1;
}
// Build up an unordered set of vertices from the VTK file. For each
// vertex a hashed value of the coordinates is generated to within a
// given tolerance.
n = vtkPoints->GetNumberOfPoints();
for (i = 0; i < n; ++i)
{
vtkPoints->GetPoint(i,p);
val = vtkPData->GetScalars(name.c_str())->GetTuple1(i);
boost::shared_ptr<Vertex> v(new Vertex(p[0],p[1],p[2],val,factor));
points.insert(v);
}
// Now process each vertex of each element in the mesh
for (i = 0; i < Exp->GetNumElmts(); ++i)
{
for (j = 0; j < e->GetNverts(); ++j)
{
// Get the index of the coefficient corresponding to this vertex
coeff_idx = Exp->GetCoeff_Offset(i) + e->GetVertexMap(j);
// Get the coordinates of the vertex
vert = e->as<LocalRegions::Expansion2D>()->GetGeom2D()
->GetVertex(j);
vert->GetCoords(x,y,z);
// Look up the vertex in the VertexSet
boost::shared_ptr<Vertex> v(new Vertex(x,y,z,0.0,factor));
vIter = points.find(v);
// If not found, maybe the tolerance should be reduced?
// If found, record the scalar value from the VTK file in the
// corresponding coefficient.
if (vIter == points.end())
{
cerr << "Vertex " << i << " not found. Looking for ("
<< x << ", " << y << ", " << z << ")" << endl;
}
else
{
Exp->UpdateCoeffs()[coeff_idx] = (*vIter)->scalar;
}
}
}
Exp->SetPhysState(false);
//-----------------------------------------------
// Write solution to file
std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
= Exp->GetFieldDefinitions();
std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
for(i = 0; i < FieldDef.size(); ++i)
{
FieldDef[i]->m_fields.push_back(outname);
Exp->AppendFieldData(FieldDef[i], FieldData[i]);
}
LibUtilities::FieldIO vFld(vSession->GetComm());
vFld.Write(outfile, FieldDef, FieldData);
//-----------------------------------------------
}
catch (...) {
cout << "An error occurred." << endl;
}
}
bool operator== ( const VertexSharedPtr v1,
const VertexSharedPtr v2 
)

Define comparison operator for the vertex struct.

Definition at line 51 of file VtkToFld.cpp.

{
return v1->operator==(*v2);
}