Nektar++
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. More...
 
- 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. More...
 
- 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. More...
 
void PrintConfig ()
 Print out all configuration options for a module. More...
 
void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
bool GetRequireEquiSpaced (void)
 
void SetRequireEquiSpaced (bool pVal)
 
void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 
 Module (MeshSharedPtr p_m)
 
void RegisterConfig (string key, string value)
 
void PrintConfig ()
 
void SetDefaults ()
 
MeshSharedPtr GetMesh ()
 
virtual void ProcessVertices ()
 Extract element vertices. More...
 

Static Public Member Functions

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

Static Public Attributes

static ModuleKey className
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::Utilities::InputModule
void PrintSummary ()
 Print summary of elements. More...
 
void PrintSummary ()
 Print summary of elements. More...
 
- Protected Member Functions inherited from Nektar::Utilities::Module
 Module ()
 
virtual void ProcessEdges (bool ReprocessEdges=true)
 Extract element edges. More...
 
virtual void ProcessFaces (bool ReprocessFaces=true)
 Extract element faces. More...
 
virtual void ProcessElements ()
 Generate element IDs. More...
 
virtual void ProcessComposites ()
 Generate composites. More...
 
void ReorderPrisms (PerMap &perFaces)
 Reorder node IDs so that prisms and tetrahedra are aligned correctly. More...
 
void PrismLines (int prism, PerMap &perFaces, set< int > &prismsDone, vector< ElementSharedPtr > &line)
 
- Protected Attributes inherited from Nektar::Utilities::InputModule
set< string > m_allowedFiles
 
std::ifstream m_mshFile
 Input stream. More...
 
- Protected Attributes inherited from Nektar::Utilities::Module
FieldSharedPtr m_f
 Field object. More...
 
map< string, ConfigOptionm_config
 List of configuration values. More...
 
bool m_requireEquiSpaced
 
MeshSharedPtr m_mesh
 Mesh object. More...
 

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.

58  : InputModule(m)
59  {
60 
61  }
Nektar::Utilities::InputVtk::~InputVtk ( )
virtual

Definition at line 63 of file InputVtk.cpp.

64  {
65 
66  }

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.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

50  {
52  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
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().

80  {
81  if (m_mesh->m_verbose)
82  {
83  cout << "InputVtk: Start reading file..." << endl;
84  }
85 
86  vtkPolyDataReader *vtkMeshReader = vtkPolyDataReader::New();
87  vtkMeshReader->SetFileName(m_config["infile"].as<string>().c_str());
88  vtkMeshReader->Update();
89  vtkPolyData *vtkMesh = vtkMeshReader->GetOutput();
90 
91  vtkPoints *vtkPoints = vtkMesh->GetPoints();
92 
93  const int numCellTypes = 3;
94  vtkCellArray* vtkCells[numCellTypes];
95  LibUtilities::ShapeType vtkCellTypes[numCellTypes];
96  int vtkNumPoints[numCellTypes];
97  vtkCells[0] = vtkMesh->GetPolys();
98  vtkCells[1] = vtkMesh->GetStrips();
99  vtkCells[2] = vtkMesh->GetLines();
100  vtkCellTypes[0] = LibUtilities::eTriangle;
101  vtkCellTypes[1] = LibUtilities::eTriangle;
102  vtkCellTypes[2] = LibUtilities::eSegment;
103  vtkNumPoints[0] = 3;
104  vtkNumPoints[1] = 3;
105  vtkNumPoints[2] = 2;
106 
107  vtkIdType npts;
108  vtkIdType *pts = 0;
109  double p[3];
110 
111  for (int i = 0; i < vtkPoints->GetNumberOfPoints(); ++i)
112  {
113  vtkPoints->GetPoint(i, p);
114 
115  if ((p[0] * p[0]) > 0.000001 && m_mesh->m_spaceDim < 1)
116  {
117  m_mesh->m_spaceDim = 1;
118  }
119  if ((p[1] * p[1]) > 0.000001 && m_mesh->m_spaceDim < 2)
120  {
121  m_mesh->m_spaceDim = 2;
122  }
123  if ((p[2] * p[2]) > 0.000001 && m_mesh->m_spaceDim < 3)
124  {
125  m_mesh->m_spaceDim = 3;
126  }
127 
128  m_mesh->m_node.push_back(boost::shared_ptr<Node>(new Node(i, p[0], p[1], p[2])));
129  }
130 
131  for (int c = 0; c < numCellTypes; ++c)
132  {
133  vtkCells[c]->InitTraversal();
134  for (int i = 0; vtkCells[c]->GetNextCell(npts, pts); ++i)
135  {
136  for (int j = 0; j < npts - vtkNumPoints[c] + 1; ++j)
137  {
138  // Create element tags
139  vector<int> tags;
140  tags.push_back(0); // composite
141  tags.push_back(vtkCellTypes[c]); // element type
142 
143  // Read element node list
144  vector<NodeSharedPtr> nodeList;
145  for (int k = j; k < j + vtkNumPoints[c]; ++k)
146  {
147  nodeList.push_back(m_mesh->m_node[pts[k]]);
148  }
149 
150  // Create element
151  ElmtConfig conf(vtkCellTypes[c],1,false,false);
153  CreateInstance(vtkCellTypes[c],
154  conf,nodeList,tags);
155 
156  // Determine mesh expansion dimension
157  if (E->GetDim() > m_mesh->m_expDim) {
158  m_mesh->m_expDim = E->GetDim();
159  }
160  m_mesh->m_element[E->GetDim()].push_back(E);
161  }
162  }
163  }
164 
165  ProcessVertices();
166  ProcessEdges();
167  ProcessFaces();
168  ProcessElements();
170  }
map< string, ConfigOption > m_config
List of configuration values.
MeshSharedPtr m_mesh
Mesh object.
virtual void ProcessEdges(bool ReprocessEdges=true)
Extract element edges.
boost::shared_ptr< Element > ElementSharedPtr
Shared pointer to an element.
Definition: MeshElements.h:63
virtual void ProcessVertices()
Extract element vertices.
virtual void ProcessElements()
Generate element IDs.
static std::string npts
Definition: InputFld.cpp:43
virtual void ProcessComposites()
Generate composites.
virtual void ProcessFaces(bool ReprocessFaces=true)
Extract element faces.
ElementFactory & GetElementFactory()

Member Data Documentation

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

Definition at line 53 of file InputVtk.h.