Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
Nektar::Utilities::InputVtk Class Reference

Converter for VTK files. More...

#include <InputVtk.h>

Inheritance diagram for Nektar::Utilities::InputVtk:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Utilities::InputVtk:
Collaboration graph
[legend]

Public Member Functions

 InputVtk (MeshSharedPtr m)
virtual ~InputVtk ()
virtual void Process ()
 Populate and validate required data structures.
- Public Member Functions inherited from Nektar::Utilities::InputModule
 InputModule (FieldSharedPtr p_m)
void AddFile (string fileType, string fileName)
 InputModule (MeshSharedPtr p_m)
void OpenStream ()
 Open a file for input.
- Public Member Functions inherited from Nektar::Utilities::Module
 Module (FieldSharedPtr p_f)
virtual void Process (po::variables_map &vm)=0
void RegisterConfig (string key, string value)
 Register a configuration option with a module.
void PrintConfig ()
 Print out all configuration options for a module.
void SetDefaults ()
 Sets default configuration options for those which have not been set.
bool GetRequireEquiSpaced (void)
void SetRequireEquiSpaced (bool pVal)
 Module (MeshSharedPtr p_m)
void RegisterConfig (string key, string value)
void PrintConfig ()
void SetDefaults ()
MeshSharedPtr GetMesh ()
virtual void ProcessVertices ()
 Extract element vertices.

Static Public Member Functions

static ModuleSharedPtr create (MeshSharedPtr m)
 Creates an instance of this class.

Static Public Attributes

static ModuleKey className

Additional Inherited Members

- Protected Member Functions inherited from Nektar::Utilities::InputModule
void PrintSummary ()
 Print summary of elements.
void PrintSummary ()
 Print summary of elements.
- Protected Attributes inherited from Nektar::Utilities::InputModule
set< string > m_allowedFiles
std::ifstream m_mshFile
 Input stream.

Detailed Description

Converter for VTK files.

Definition at line 46 of file InputVtk.h.

Constructor & Destructor Documentation

Nektar::Utilities::InputVtk::InputVtk ( MeshSharedPtr  m)

Definition at line 58 of file InputVtk.cpp.

{
}
Nektar::Utilities::InputVtk::~InputVtk ( )
virtual

Definition at line 63 of file InputVtk.cpp.

{
}

Member Function Documentation

static ModuleSharedPtr Nektar::Utilities::InputVtk::create ( MeshSharedPtr  m)
inlinestatic

Creates an instance of this class.

Definition at line 50 of file InputVtk.h.

{
return MemoryManager<InputVtk>::AllocateSharedPtr(m);
}
void Nektar::Utilities::InputVtk::Process ( )
virtual

Populate and validate required data structures.

Gmsh file contains a list of nodes and their coordinates, along with a list of elements and those nodes which define them. We read in and store the list of nodes in #m_node and store the list of elements in #m_element. Each new element is supplied with a list of entries from #m_node which defines the element. Finally some mesh statistics are printed.

Parameters
pFilenameFilename of Gmsh file to read.

Implements Nektar::Utilities::Module.

Definition at line 79 of file InputVtk.cpp.

References Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTriangle, Nektar::Utilities::GetElementFactory(), Nektar::Utilities::Module::m_config, Nektar::Utilities::Module::m_mesh, npts, Nektar::Utilities::Module::ProcessComposites(), Nektar::Utilities::Module::ProcessEdges(), Nektar::Utilities::Module::ProcessElements(), Nektar::Utilities::Module::ProcessFaces(), and Nektar::Utilities::Module::ProcessVertices().

{
if (m_mesh->m_verbose)
{
cout << "InputVtk: Start reading file..." << endl;
}
vtkPolyDataReader *vtkMeshReader = vtkPolyDataReader::New();
vtkMeshReader->SetFileName(m_config["infile"].as<string>().c_str());
vtkMeshReader->Update();
vtkPolyData *vtkMesh = vtkMeshReader->GetOutput();
vtkPoints *vtkPoints = vtkMesh->GetPoints();
const int numCellTypes = 3;
vtkCellArray* vtkCells[numCellTypes];
LibUtilities::ShapeType vtkCellTypes[numCellTypes];
int vtkNumPoints[numCellTypes];
vtkCells[0] = vtkMesh->GetPolys();
vtkCells[1] = vtkMesh->GetStrips();
vtkCells[2] = vtkMesh->GetLines();
vtkCellTypes[0] = LibUtilities::eTriangle;
vtkCellTypes[1] = LibUtilities::eTriangle;
vtkCellTypes[2] = LibUtilities::eSegment;
vtkNumPoints[0] = 3;
vtkNumPoints[1] = 3;
vtkNumPoints[2] = 2;
vtkIdType npts;
vtkIdType *pts = 0;
double p[3];
for (int i = 0; i < vtkPoints->GetNumberOfPoints(); ++i)
{
vtkPoints->GetPoint(i, p);
if ((p[0] * p[0]) > 0.000001 && m_mesh->m_spaceDim < 1)
{
m_mesh->m_spaceDim = 1;
}
if ((p[1] * p[1]) > 0.000001 && m_mesh->m_spaceDim < 2)
{
m_mesh->m_spaceDim = 2;
}
if ((p[2] * p[2]) > 0.000001 && m_mesh->m_spaceDim < 3)
{
m_mesh->m_spaceDim = 3;
}
m_mesh->m_node.push_back(boost::shared_ptr<Node>(new Node(i, p[0], p[1], p[2])));
}
for (int c = 0; c < numCellTypes; ++c)
{
vtkCells[c]->InitTraversal();
for (int i = 0; vtkCells[c]->GetNextCell(npts, pts); ++i)
{
for (int j = 0; j < npts - vtkNumPoints[c] + 1; ++j)
{
// Create element tags
vector<int> tags;
tags.push_back(0); // composite
tags.push_back(vtkCellTypes[c]); // element type
// Read element node list
vector<NodeSharedPtr> nodeList;
for (int k = j; k < j + vtkNumPoints[c]; ++k)
{
nodeList.push_back(m_mesh->m_node[pts[k]]);
}
// Create element
ElmtConfig conf(vtkCellTypes[c],1,false,false);
CreateInstance(vtkCellTypes[c],
conf,nodeList,tags);
// Determine mesh expansion dimension
if (E->GetDim() > m_mesh->m_expDim) {
m_mesh->m_expDim = E->GetDim();
}
m_mesh->m_element[E->GetDim()].push_back(E);
}
}
}
}

Member Data Documentation

ModuleKey Nektar::Utilities::InputVtk::className
static
Initial value:

Definition at line 53 of file InputVtk.h.