Nektar++
utilities/NekMesh/OutputModules/OutputVtk.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: OutputVtk.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: VTK file format output.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 
37 
38 #include <vtkPolyDataWriter.h>
39 #include <vtkPolyData.h>
40 #include <vtkPoints.h>
41 #include <vtkCellArray.h>
42 
43 #include "OutputVtk.h"
44 
45 using namespace std;
46 using namespace Nektar::NekMeshUtils;
47 
48 namespace Nektar
49 {
50 namespace Utilities
51 {
52 
53 ModuleKey OutputVtk::className = GetModuleFactory().RegisterCreatorFunction(
54  ModuleKey(eOutputModule, "vtk"), OutputVtk::create, "Writes a VTK file.");
55 
56 OutputVtk::OutputVtk(MeshSharedPtr m) : OutputModule(m)
57 {
58 }
59 
61 {
62 }
63 
65 {
66  if (m_mesh->m_verbose)
67  {
68  cout << "OutputVtk: Writing file..." << endl;
69  }
70 
71  vtkPolyData *vtkMesh = vtkPolyData::New();
72  vtkPoints *vtkPoints = vtkPoints::New();
73  vtkCellArray *vtkPolys = vtkCellArray::New();
74 
75  std::set<NodeSharedPtr> tmp(m_mesh->m_vertexSet.begin(),
76  m_mesh->m_vertexSet.end());
77 
78  for (auto &n : tmp)
79  {
80  vtkPoints->InsertPoint(n->m_id, n->m_x, n->m_y, n->m_z);
81  }
82 
83  vtkIdType p[8];
84  vector<ElementSharedPtr> &elmt = m_mesh->m_element[m_mesh->m_expDim];
85  for (int i = 0; i < elmt.size(); ++i)
86  {
87  int vertexCount = elmt[i]->GetVertexCount();
88  for (int j = 0; j < vertexCount; ++j)
89  {
90  p[j] = elmt[i]->GetVertex(j)->m_id;
91  }
92  vtkPolys->InsertNextCell(vertexCount, &p[0]);
93  }
94 
95  vtkMesh->SetPoints(vtkPoints);
96  vtkMesh->SetPolys(vtkPolys);
97 
98  // Write out the new mesh
99  vtkPolyDataWriter *vtkMeshWriter = vtkPolyDataWriter::New();
100  vtkMeshWriter->SetFileName(m_config["outfile"].as<string>().c_str());
101 #if VTK_MAJOR_VERSION <= 5
102  vtkMeshWriter->SetInput(vtkMesh);
103 #else
104  vtkMeshWriter->SetInputData(vtkMesh);
105 #endif
106  vtkMeshWriter->Update();
107 }
108 }
109 }
Abstract base class for output modules.
STL namespace.
std::shared_ptr< Mesh > MeshSharedPtr
Shared pointer to a mesh.
Definition: Mesh.h:156
std::pair< ModuleType, std::string > ModuleKey
virtual void Process()
Write mesh to output file.
std::map< std::string, ConfigOption > m_config
List of configuration values.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
std::pair< ModuleType, std::string > ModuleKey
ModuleFactory & GetModuleFactory()