Chaste  Release::2017.1
VtkMeshReader.cpp
1 /*
2 
3 Copyright (C) Fujitsu Laboratories of Europe, 2009
4 
5 */
6 
7 /*
8 
9 Copyright (c) 2005-2017, University of Oxford.
10 All rights reserved.
11 
12 University of Oxford means the Chancellor, Masters and Scholars of the
13 University of Oxford, having an administrative office at Wellington
14 Square, Oxford OX1 2JD, UK.
15 
16 This file is part of Chaste.
17 
18 Redistribution and use in source and binary forms, with or without
19 modification, are permitted provided that the following conditions are met:
20  * Redistributions of source code must retain the above copyright notice,
21  this list of conditions and the following disclaimer.
22  * Redistributions in binary form must reproduce the above copyright notice,
23  this list of conditions and the following disclaimer in the documentation
24  and/or other materials provided with the distribution.
25  * Neither the name of the University of Oxford nor the names of its
26  contributors may be used to endorse or promote products derived from this
27  software without specific prior written permission.
28 
29 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
30 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
31 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
32 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
33 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
34 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
35 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
36 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
37 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
38 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
39 
40 */
41 
42 #ifdef CHASTE_VTK
43 
44 #include "VtkMeshReader.hpp"
45 #include "Exception.hpp"
46 
47 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
49  mIndexFromZero(true),
50  mNumNodes(0),
51  mNumElements(0),
52  mNumFaces(0),
53  mNumCableElements(0),
54  mNodesRead(0),
55  mElementsRead(0),
56  mFacesRead(0),
57  mBoundaryFacesRead(0),
58  mBoundaryFacesSkipped(0),
59  mCableElementsRead(0),
60  mNumElementAttributes(0),
61  mNumFaceAttributes(0),
62  mNumCableElementAttributes(0),
63  mOrderOfElements(1),
64  mNodesPerElement(4),
65  mVtkCellType(VTK_TETRA)
66 {
67  vtkXMLUnstructuredGridReader* vtk_xml_unstructured_grid_reader;
68 
69  // Check file exists
70  mVtuFile.open(pathBaseName.c_str());
71  if (!mVtuFile.is_open())
72  {
73  EXCEPTION("Could not open VTU file: " + pathBaseName);
74  }
75  mVtuFile.close();
76 
77  // Load the mesh geometry and data from a file
78  vtk_xml_unstructured_grid_reader = vtkXMLUnstructuredGridReader::New();
79  vtk_xml_unstructured_grid_reader->SetFileName( pathBaseName.c_str() );
80  vtk_xml_unstructured_grid_reader->Update();
81  mpVtkUnstructuredGrid = vtk_xml_unstructured_grid_reader->GetOutput();
83  vtk_xml_unstructured_grid_reader->Delete();
84 }
85 
86 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
88 {
89 
90  mNumNodes = mpVtkUnstructuredGrid->GetNumberOfPoints();
91  unsigned num_cells = mpVtkUnstructuredGrid->GetNumberOfCells();
92 
93  if (ELEMENT_DIM == 2u)
94  {
95  mNodesPerElement = 3;
96  mVtkCellType = VTK_TRIANGLE;
97  }
98  else if (ELEMENT_DIM == 1u)
99  {
100  mNodesPerElement = 2;
101  mVtkCellType = VTK_LINE;
102  }
103 
104  //Determine if we have multiple cell types - such as cable elements in addition to tets/triangles
105  vtkCellTypes* cell_types = vtkCellTypes::New();
106  mpVtkUnstructuredGrid->GetCellTypes(cell_types);
107 
108  if (cell_types->GetNumberOfTypes() > 1)
109  {
111  for (unsigned cell_id = 0; cell_id < num_cells; ++cell_id)
112  {
113  if (mpVtkUnstructuredGrid->GetCellType(cell_id) == mVtkCellType)
114  {
115  ++mNumElements;
116  assert(mNumCableElements == 0); //We expect all the simplices first and then the cables at the end of the array
117  }
118  else if (mpVtkUnstructuredGrid->GetCellType(cell_id) == VTK_LINE)
119  {
121  }
122  else
123  {
125  }
126  }
127  }
128  else
129  {
130  //There is only 1 cell type, so all cells are elements
131  mNumElements = num_cells;
132  }
133 
134  cell_types->Delete();
135 
136 
137  // Extract the surface faces
138  if (ELEMENT_DIM == 2u)
139  {
140  vtkDataSetSurfaceFilter* p_surface = vtkDataSetSurfaceFilter::New();
141  mpVtkFilterEdges = vtkFeatureEdges::New();
142 #if VTK_MAJOR_VERSION >= 6
143  p_surface->SetInputData(mpVtkUnstructuredGrid);
144  mpVtkFilterEdges->SetInputConnection(p_surface->GetOutputPort());
145 #else
146  p_surface->SetInput(mpVtkUnstructuredGrid);
147  mpVtkFilterEdges->SetInput(p_surface->GetOutput());
148 #endif
149  mpVtkFilterEdges->Update();
150  mNumFaces = mpVtkFilterEdges->GetOutput()->GetNumberOfCells();
151  p_surface->Delete();
152  }
153  else if (ELEMENT_DIM == 3u)
154  {
155  mpVtkGeometryFilter = vtkGeometryFilter::New();
156 #if VTK_MAJOR_VERSION >= 6
158 #else
160 #endif
161  mpVtkGeometryFilter->Update();
162 
163  mNumFaces = mpVtkGeometryFilter->GetOutput()->GetNumberOfCells();
164  if (mNumCableElements > 0)
165  {
166  //The boundary face filter includes the cable elements - get rid of them
167  unsigned num_all_cells = mNumFaces;
168  for (unsigned i=0; i<num_all_cells; i++)
169  {
170  if (mpVtkGeometryFilter->GetOutput()->GetCellType(i) == VTK_LINE)
171  {
172  mNumFaces--;
173  }
174  }
175  }
176  }
177  else if (ELEMENT_DIM == 1u)
178  {
179  mNumFaces = 0;
180  vtkSmartPointer<vtkIdList> enclosing_cells = vtkSmartPointer<vtkIdList>::New();
181  assert(mNumNodes == (unsigned) mpVtkUnstructuredGrid->GetNumberOfPoints());
182  for (unsigned point_index = 0; point_index < mNumNodes; ++point_index)
183  {
184  mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
185 
186  if (enclosing_cells->GetNumberOfIds() == 1u)
187  {
188  mNumFaces++;
189  }
190  }
191  }
192  else
193  {
195  }
196 }
197 
198 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
199 VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::VtkMeshReader(vtkUnstructuredGrid* p_vtkUnstructuredGrid) :
200  mIndexFromZero(true),
201  mNumNodes(0),
202  mNumElements(0),
203  mNumFaces(0),
205  mNodesRead(0),
206  mElementsRead(0),
207  mFacesRead(0),
214  mOrderOfElements(1),
215  mNodesPerElement(4),
216  mVtkCellType(VTK_TETRA)
217 {
218  mpVtkUnstructuredGrid = p_vtkUnstructuredGrid;
220 }
221 
222 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
224 {
225  if (ELEMENT_DIM == 3u)
226  {
227  mpVtkGeometryFilter->Delete();
228  }
229  else if (ELEMENT_DIM == 2u)
230  {
231  mpVtkFilterEdges->Delete();
232  }
233 }
234 
235 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
237 {
238  return mNumElements;
239 }
240 
241 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
243 {
244  return mNumCableElements;
245 }
246 
247 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
249 {
250  return mNumNodes;
251 }
252 
253 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
255 {
256  return mNumFaces;
257 }
258 
259 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
261 {
262  return mNumFaces;
263 }
264 
265 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
267 {
268  return mNumElementAttributes;
269 }
270 
271 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
273 {
275 }
276 
277 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
279 {
280  return mNumFaceAttributes;
281 }
282 
283 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
285 {
286  mNodesRead = 0;
287  mElementsRead = 0;
288  mFacesRead = 0;
289  mBoundaryFacesRead = 0;
291  mCableElementsRead = 0;
292 }
293 
294 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
296 {
297  if (mNodesRead >= mNumNodes)
298  {
299  EXCEPTION( "Trying to read data for a node that doesn't exist" );
300  }
301 
302  std::vector<double> next_node;
303 
304  for (unsigned i = 0; i < 3; i++)
305  {
306  next_node.push_back( mpVtkUnstructuredGrid->GetPoint(mNodesRead)[i] );
307  }
308 
309  mNodesRead++;
310  return next_node;
311 }
312 
313 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
315 {
317  {
318  EXCEPTION( "Trying to read data for an element that doesn't exist" );
319  }
320 
322  {
323  EXCEPTION("Element is not of expected type (vtkTetra/vtkTriangle)");
324  }
325 
326  ElementData next_element_data;
327 
328  for (unsigned i = 0; i < mNodesPerElement; i++)
329  {
330  next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(mElementsRead)->GetPointId(i));
331  }
332 
333  // \todo implement method to read element data properly (currently returns zero always...)
334  next_element_data.AttributeValue = 0;
335 
336  mElementsRead++;
337  return next_element_data;
338 }
339 
340 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
342 {
344  {
345  EXCEPTION( "Trying to read data for a cable element that doesn't exist" );
346  }
347 
348  unsigned next_index = mNumElements + mCableElementsRead;
349  assert(mpVtkUnstructuredGrid->GetCellType(next_index)==VTK_LINE);
350 
351  ElementData next_element_data;
352 
353  for (unsigned i = 0; i < 2; i++)
354  {
355  next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(next_index)->GetPointId(i));
356  }
357 
358  vtkDataArray *p_scalars = mpVtkUnstructuredGrid->GetCellData()->GetArray( "Cable radius" );
359  next_element_data.AttributeValue = p_scalars->GetTuple(next_index)[0];
360 
361  mCableElementsRead++;
362  return next_element_data;
363 }
364 
365 
366 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
368 {
370  {
371  EXCEPTION( "Trying to read data for a boundary element that doesn't exist");
372  }
373 
374  ElementData next_face_data;
375 
376  if (ELEMENT_DIM == 3u)
377  {
378  while (mpVtkGeometryFilter->GetOutput()->GetCellType(mBoundaryFacesRead + mBoundaryFacesSkipped) == VTK_LINE)
379  {
381  }
382  for (unsigned i = 0; i < (mNodesPerElement-1); i++)
383  {
384  next_face_data.NodeIndices.push_back(mpVtkGeometryFilter->GetOutput()->GetCell(mBoundaryFacesRead + mBoundaryFacesSkipped)->GetPointId(i));
385  }
386  }
387  else if (ELEMENT_DIM == 2u)
388  {
389  for (unsigned i = 0; i < (mNodesPerElement-1); i++)
390  {
391  next_face_data.NodeIndices.push_back(mpVtkFilterEdges->GetOutput()->GetCell(mBoundaryFacesRead)->GetPointId(i));
392  }
393  }
394  else if (ELEMENT_DIM == 1u)
395  {
396  vtkSmartPointer<vtkIdList> enclosing_cells = vtkSmartPointer<vtkIdList>::New();
397  assert(mNumNodes == (unsigned) mpVtkUnstructuredGrid->GetNumberOfPoints());
398  for (unsigned point_index = mBoundaryFacesSkipped; point_index < mNumNodes; ++point_index)
399  {
401  mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
402 
403  if (enclosing_cells->GetNumberOfIds() == 1u)
404  {
405  next_face_data.NodeIndices.push_back(point_index);
406  break;
407  }
408  }
409  }
410  else
411  {
413  }
414 
416  return next_face_data;
417 }
418 
419 
420 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
421 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<double>& dataPayload)
422 {
423  vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
424 
425  if (!p_cell_data->HasArray(dataName.c_str()))
426  {
427  EXCEPTION("No cell data '" + dataName + "'");
428  }
429 
430  vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
431  if (p_scalars->GetNumberOfComponents() != 1)
432  {
433  EXCEPTION("The cell data '" + dataName + "' is not scalar data.");
434  }
435 
436  dataPayload.clear();
437  for (unsigned i = 0; i < mNumElements; i++)
438  {
439  dataPayload.push_back( p_scalars->GetTuple(i)[0] );
440  }
441 }
442 
443 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
444 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
445 {
446  vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
447 
448  if (!p_cell_data->HasArray(dataName.c_str()))
449  {
450  EXCEPTION("No cell data '" + dataName + "'");
451  }
452 
453  vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
454  if (p_scalars->GetNumberOfComponents() != 3)
455  {
456  EXCEPTION("The cell data '" + dataName + "' is not 3-vector data.");
457  }
458 
459  dataPayload.clear();
460  for (unsigned i = 0; i < mNumElements; i++)
461  {
462  c_vector <double, SPACE_DIM> data;
463  for (unsigned j = 0; j < SPACE_DIM; j++)
464  {
465  data[j]= p_scalars->GetTuple(i)[j];
466  }
467  dataPayload.push_back(data);
468  }
469 }
470 
471 
472 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
473 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<double>& dataPayload)
474 {
475  vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
476 
477  if (!p_point_data->HasArray(dataName.c_str()))
478  {
479  EXCEPTION("No point data '" + dataName + "'");
480  }
481 
482  vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
483  if (p_scalars->GetNumberOfComponents() != 1)
484  {
485  EXCEPTION("The point data '" + dataName + "' is not scalar data.");
486  }
487 
488  dataPayload.clear();
489 
490  for (unsigned i = 0; i < mNumNodes; i++)
491  {
492  dataPayload.push_back( p_scalars->GetTuple(i)[0] );
493  }
494 }
495 
496 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
497 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
498 {
499  vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
500 
501  if (!p_point_data->HasArray(dataName.c_str()))
502  {
503  EXCEPTION("No point data '" + dataName + "'");
504  }
505 
506  vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
507 
508  if (p_scalars->GetNumberOfComponents() != 3)
509  {
510  EXCEPTION("The point data '" + dataName + "' is not 3-vector data.");
511  }
512  dataPayload.clear();
513  for (unsigned i = 0; i < mNumNodes; i++)
514  {
515  c_vector<double, SPACE_DIM> data;
516  for (unsigned j = 0; j< SPACE_DIM; j++)
517  {
518  data[j] = p_scalars->GetTuple(i)[j];
519  }
520  dataPayload.push_back( data );
521  }
522 }
523 
524 
525 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
527 {
528  return mpVtkUnstructuredGrid;
529 }
530 
531 // Explicit instantiation
532 template class VtkMeshReader<1,1>;
533 template class VtkMeshReader<1,2>;
534 template class VtkMeshReader<1,3>;
535 template class VtkMeshReader<2,2>;
536 template class VtkMeshReader<2,3>;
537 template class VtkMeshReader<3,3>;
538 #endif // CHASTE_VTK
void CommonConstructor()
unsigned GetNumCableElementAttributes() const
std::ifstream mVtuFile
unsigned mBoundaryFacesSkipped
unsigned mNodesRead
unsigned mNumNodes
#define EXCEPTION(message)
Definition: Exception.hpp:143
unsigned GetNumEdges() const
unsigned mNumFaceAttributes
void GetPointData(std::string dataName, std::vector< double > &dataPayload)
unsigned GetNumFaceAttributes() const
unsigned mOrderOfElements
vtkGeometryFilter * mpVtkGeometryFilter
ElementData GetNextFaceData()
unsigned GetNumElementAttributes() const
#define NEVER_REACHED
Definition: Exception.hpp:206
ElementData GetNextElementData()
unsigned mNumElementAttributes
std::vector< unsigned > NodeIndices
unsigned mNumCableElements
unsigned GetNumElements() const
virtual ~VtkMeshReader()
unsigned mNodesPerElement
unsigned mNumFaces
ElementData GetNextCableElementData()
unsigned mCableElementsRead
unsigned GetNumCableElements() const
unsigned mBoundaryFacesRead
vtkFeatureEdges * mpVtkFilterEdges
void GetCellData(std::string dataName, std::vector< double > &dataPayload)
unsigned mNumElements
unsigned mNumCableElementAttributes
unsigned GetNumFaces() const
unsigned GetNumNodes() const
vtkUnstructuredGrid * OutputMeshAsVtkUnstructuredGrid()
unsigned mFacesRead
std::vector< double > GetNextNode()
vtkSmartPointer< vtkUnstructuredGrid > mpVtkUnstructuredGrid
VtkMeshReader(std::string pathBaseName)
unsigned mElementsRead