Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
InputVtk.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: InputVtk.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 converter.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <string>
37 #include <fstream>
38 #include <iostream>
39 using namespace std;
40 
41 #include <vtkPolyDataReader.h>
42 #include <vtkPolyData.h>
43 #include <vtkPoints.h>
44 #include <vtkCellArray.h>
45 
46 #include "MeshElements.h"
47 #include "InputVtk.h"
48 
49 namespace Nektar
50 {
51  namespace Utilities
52  {
53  ModuleKey InputVtk::className =
55  ModuleKey(eInputModule, "vtk"), InputVtk::create,
56  "Reads VTK format.");
57 
58  InputVtk::InputVtk(MeshSharedPtr m) : InputModule(m)
59  {
60 
61  }
62 
64  {
65 
66  }
67 
68 
69  /**
70  * Gmsh file contains a list of nodes and their coordinates, along with
71  * a list of elements and those nodes which define them. We read in and
72  * store the list of nodes in #m_node and store the list of elements in
73  * #m_element. Each new element is supplied with a list of entries from
74  * #m_node which defines the element. Finally some mesh statistics are
75  * printed.
76  *
77  * @param pFilename Filename of Gmsh file to read.
78  */
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  }
171  }
172 }