Nektar++
MeshConvert/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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: VTK file format output.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <set>
37 #include <string>
38 using namespace std;
39 
40 #include <vtkPolyDataWriter.h>
41 #include <vtkPolyData.h>
42 #include <vtkPoints.h>
43 #include <vtkCellArray.h>
44 
45 #include "../MeshElements.h"
46 #include "OutputVtk.h"
47 
48 namespace Nektar
49 {
50  namespace Utilities
51  {
52  ModuleKey OutputVtk::className =
54  ModuleKey(eOutputModule, "vtk"), OutputVtk::create,
55  "Writes a VTK file.");
56 
57  OutputVtk::OutputVtk(MeshSharedPtr m) : OutputModule(m)
58  {
59 
60  }
61 
63  {
64 
65  }
66 
68  {
69  if (m_mesh->m_verbose)
70  {
71  cout << "OutputVtk: Writing file..." << endl;
72  }
73 
74  vtkPolyData *vtkMesh = vtkPolyData::New();
75  vtkPoints *vtkPoints = vtkPoints::New();
76  vtkCellArray *vtkPolys = vtkCellArray::New();
77 
79 
80  std::set<NodeSharedPtr> tmp(
81  m_mesh->m_vertexSet.begin(),
82  m_mesh->m_vertexSet.end());
83 
84  for (it = tmp.begin(); it != tmp.end(); ++it)
85  {
86  NodeSharedPtr n = *it;
87  vtkPoints->InsertPoint(n->m_id, n->m_x, n->m_y, n->m_z);
88  }
89 
90  vtkIdType p[8];
91  vector<ElementSharedPtr> &elmt =
92  m_mesh->m_element[m_mesh->m_expDim];
93  for(int i = 0; i < elmt.size(); ++i)
94  {
95  int vertexCount = elmt[i]->GetVertexCount();
96  for (int j = 0; j < vertexCount; ++j)
97  {
98  p[j] = elmt[i]->GetVertex(j)->m_id;
99  }
100  vtkPolys->InsertNextCell(vertexCount, &p[0]);
101  }
102 
103  vtkMesh->SetPoints(vtkPoints);
104  vtkMesh->SetPolys(vtkPolys);
105 
106  // Write out the new mesh
107  vtkPolyDataWriter *vtkMeshWriter = vtkPolyDataWriter::New();
108  vtkMeshWriter->SetFileName(m_config["outfile"].as<string>().c_str());
109 #if VTK_MAJOR_VERSION <= 5
110  vtkMeshWriter->SetInput(vtkMesh);
111 #else
112  vtkMeshWriter->SetInputData(vtkMesh);
113 #endif
114  vtkMeshWriter->Update();
115  }
116  }
117 }
pair< ModuleType, string > ModuleKey
Abstract base class for output modules.
map< string, ConfigOption > m_config
List of configuration values.
STL namespace.
MeshSharedPtr m_mesh
Mesh object.
boost::shared_ptr< Node > NodeSharedPtr
Shared pointer to a Node.
Definition: MeshElements.h:195
virtual void Process()
Write mesh to output file.
boost::shared_ptr< Mesh > MeshSharedPtr
Shared pointer to a mesh.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
ModuleFactory & GetModuleFactory()
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215