Chaste  Release::3.4
VtkMeshReader.cpp
1 /*
2 
3 Copyright (C) Fujitsu Laboratories of Europe, 2009
4 
5 */
6 
7 /*
8 
9 Copyright (c) 2005-2016, 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 
43 
44 #ifdef CHASTE_VTK
45 
46 #include "VtkMeshReader.hpp"
47 #include "Exception.hpp"
48 
49 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
51  mIndexFromZero(true),
52  mNumNodes(0),
53  mNumElements(0),
54  mNumFaces(0),
55  mNumCableElements(0),
56  mNodesRead(0),
57  mElementsRead(0),
58  mFacesRead(0),
59  mBoundaryFacesRead(0),
60  mBoundaryFacesSkipped(0),
61  mCableElementsRead(0),
62  mNumElementAttributes(0),
63  mNumFaceAttributes(0),
64  mNumCableElementAttributes(0),
65  mOrderOfElements(1),
66  mNodesPerElement(4),
67  mVtkCellType(VTK_TETRA)
68 {
69  vtkXMLUnstructuredGridReader* vtk_xml_unstructured_grid_reader;
70 
71  // Check file exists
72  mVtuFile.open(pathBaseName.c_str());
73  if ( !mVtuFile.is_open() )
74  {
75  EXCEPTION("Could not open VTU file: " + pathBaseName);
76  }
77  mVtuFile.close();
78 
79  // Load the mesh geometry and data from a file
80  vtk_xml_unstructured_grid_reader = vtkXMLUnstructuredGridReader::New();
81  vtk_xml_unstructured_grid_reader->SetFileName( pathBaseName.c_str() );
82  vtk_xml_unstructured_grid_reader->Update();
83  mpVtkUnstructuredGrid = vtk_xml_unstructured_grid_reader->GetOutput();
85  vtk_xml_unstructured_grid_reader->Delete();
86 }
87 
88 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
90 {
91 
92  mNumNodes = mpVtkUnstructuredGrid->GetNumberOfPoints();
93  unsigned num_cells = mpVtkUnstructuredGrid->GetNumberOfCells();
94 
95  if (ELEMENT_DIM == 2u)
96  {
97  mNodesPerElement = 3;
98  mVtkCellType = VTK_TRIANGLE;
99  }
100  else if (ELEMENT_DIM == 1u)
101  {
102  mNodesPerElement = 2;
103  mVtkCellType = VTK_LINE;
104  }
105 
106  //Determine if we have multiple cell types - such as cable elements in addition to tets/triangles
107  vtkCellTypes* cell_types = vtkCellTypes::New();
108  mpVtkUnstructuredGrid->GetCellTypes(cell_types);
109 
110  if(cell_types->GetNumberOfTypes() > 1)
111  {
112  mNumCableElementAttributes = 1;
113  for(unsigned cell_id = 0; cell_id < num_cells; ++cell_id)
114  {
115  if(mpVtkUnstructuredGrid->GetCellType(cell_id) == mVtkCellType)
116  {
117  ++mNumElements;
118  assert(mNumCableElements == 0); //We expect all the simplices first and then the cables at the end of the array
119  }
120  else if(mpVtkUnstructuredGrid->GetCellType(cell_id) == VTK_LINE)
121  {
122  ++mNumCableElements;
123  }
124  else
125  {
127  }
128  }
129  }
130  else
131  {
132  //There is only 1 cell type, so all cells are elements
133  mNumElements = num_cells;
134  }
135 
136  cell_types->Delete();
137 
138 
139  // Extract the surface faces
140  if (ELEMENT_DIM == 2u)
141  {
142  vtkDataSetSurfaceFilter* p_surface = vtkDataSetSurfaceFilter::New();
143  mpVtkFilterEdges = vtkFeatureEdges::New();
144 #if VTK_MAJOR_VERSION >= 6
145  p_surface->SetInputData(mpVtkUnstructuredGrid);
146  mpVtkFilterEdges->SetInputConnection(p_surface->GetOutputPort());
147 #else
148  p_surface->SetInput(mpVtkUnstructuredGrid);
149  mpVtkFilterEdges->SetInput(p_surface->GetOutput());
150 #endif
151  mpVtkFilterEdges->Update();
152  mNumFaces = mpVtkFilterEdges->GetOutput()->GetNumberOfCells();
153  p_surface->Delete();
154  }
155  else if (ELEMENT_DIM == 3u)
156  {
157  mpVtkGeometryFilter = vtkGeometryFilter::New();
158 #if VTK_MAJOR_VERSION >= 6
159  mpVtkGeometryFilter->SetInputData(mpVtkUnstructuredGrid);
160 #else
161  mpVtkGeometryFilter->SetInput(mpVtkUnstructuredGrid);
162 #endif
163  mpVtkGeometryFilter->Update();
164 
165  mNumFaces = mpVtkGeometryFilter->GetOutput()->GetNumberOfCells();
166  if (mNumCableElements > 0)
167  {
168  //The boundary face filter includes the cable elements - get rid of them
169  unsigned num_all_cells = mNumFaces;
170  for (unsigned i=0; i<num_all_cells; i++)
171  {
172  if (mpVtkGeometryFilter->GetOutput()->GetCellType(i) == VTK_LINE)
173  {
174  mNumFaces--;
175  }
176  }
177  }
178  }
179  else if (ELEMENT_DIM == 1u)
180  {
181  mNumFaces = 0;
182  vtkSmartPointer<vtkIdList> enclosing_cells = vtkSmartPointer<vtkIdList>::New();
183  assert(mNumNodes == (unsigned) mpVtkUnstructuredGrid->GetNumberOfPoints());
184  for (unsigned point_index = 0; point_index < mNumNodes; ++point_index)
185  {
186  mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
187 
188  if (enclosing_cells->GetNumberOfIds() == 1u)
189  {
190  mNumFaces++;
191  }
192  }
193  }
194  else
195  {
197  }
198 }
199 
200 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
201 VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::VtkMeshReader(vtkUnstructuredGrid* p_vtkUnstructuredGrid) :
202  mIndexFromZero(true),
203  mNumNodes(0),
204  mNumElements(0),
205  mNumFaces(0),
206  mNumCableElements(0),
207  mNodesRead(0),
208  mElementsRead(0),
209  mFacesRead(0),
210  mBoundaryFacesRead(0),
211  mBoundaryFacesSkipped(0),
212  mCableElementsRead(0),
213  mNumElementAttributes(0),
214  mNumFaceAttributes(0),
215  mNumCableElementAttributes(0),
216  mOrderOfElements(1),
217  mNodesPerElement(4),
218  mVtkCellType(VTK_TETRA)
219 {
220  mpVtkUnstructuredGrid = p_vtkUnstructuredGrid;
222 }
223 
224 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
226 {
227  if (ELEMENT_DIM == 3u)
228  {
229  mpVtkGeometryFilter->Delete();
230  }
231  else if (ELEMENT_DIM == 2u)
232  {
233  mpVtkFilterEdges->Delete();
234  }
235 }
236 
237 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
239 {
240  return mNumElements;
241 }
242 
243 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
245 {
246  return mNumCableElements;
247 }
248 
249 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
251 {
252  return mNumNodes;
253 }
254 
255 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
257 {
258  return mNumFaces;
259 }
260 
261 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
263 {
264  return mNumFaces;
265 }
266 
267 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
269 {
270  return mNumElementAttributes;
271 }
272 
273 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
275 {
276  return mNumCableElementAttributes;
277 }
278 
279 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
281 {
282  return mNumFaceAttributes;
283 }
284 
285 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
287 {
288  mNodesRead=0;
289  mElementsRead=0;
290  mFacesRead=0;
291  mBoundaryFacesRead=0;
292  mBoundaryFacesSkipped=0;
293  mCableElementsRead=0;
294 }
295 
296 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
298 {
299  mpVtkUnstructuredGrid->Initialize();
300 }
301 
302 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
304 {
305  if ( mNodesRead >= mNumNodes )
306  {
307  EXCEPTION( "Trying to read data for a node that doesn't exist" );
308  }
309 
310  std::vector<double> next_node;
311 
312  for (unsigned i = 0; i < 3; i++)
313  {
314  next_node.push_back( mpVtkUnstructuredGrid->GetPoint(mNodesRead)[i] );
315  }
316 
317  mNodesRead++;
318  return next_node;
319 }
320 
321 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
323 {
324  if ( mElementsRead >= mNumElements )
325  {
326  EXCEPTION( "Trying to read data for an element that doesn't exist" );
327  }
328 
329  if (mpVtkUnstructuredGrid->GetCellType(mElementsRead)!= mVtkCellType)
330  {
331  EXCEPTION("Element is not of expected type (vtkTetra/vtkTriangle)");
332  }
333 
334  ElementData next_element_data;
335 
336  for (unsigned i = 0; i < mNodesPerElement; i++)
337  {
338  next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(mElementsRead)->GetPointId(i));
339  }
340 
341  // \todo implement method to read element data properly (currently returns zero always...)
342  next_element_data.AttributeValue = 0;
343 
344  mElementsRead++;
345  return next_element_data;
346 }
347 
348 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
350 {
351  if ( mCableElementsRead >= mNumCableElements )
352  {
353  EXCEPTION( "Trying to read data for a cable element that doesn't exist" );
354  }
355 
356  unsigned next_index = mNumElements + mCableElementsRead;
357  assert(mpVtkUnstructuredGrid->GetCellType(next_index)==VTK_LINE);
358 
359  ElementData next_element_data;
360 
361  for (unsigned i = 0; i < 2; i++)
362  {
363  next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(next_index)->GetPointId(i));
364  }
365 
366  vtkDataArray *p_scalars = mpVtkUnstructuredGrid->GetCellData()->GetArray( "Cable radius" );
367  next_element_data.AttributeValue = p_scalars->GetTuple(next_index)[0];
368 
369  mCableElementsRead++;
370  return next_element_data;
371 }
372 
373 
374 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
376 {
377  if ( mBoundaryFacesRead >= mNumFaces)
378  {
379  EXCEPTION( "Trying to read data for a boundary element that doesn't exist");
380  }
381 
382  ElementData next_face_data;
383 
384  if (ELEMENT_DIM == 3u)
385  {
386  while (mpVtkGeometryFilter->GetOutput()->GetCellType(mBoundaryFacesRead + mBoundaryFacesSkipped) == VTK_LINE)
387  {
388  mBoundaryFacesSkipped++;
389  }
390  for (unsigned i = 0; i < (mNodesPerElement-1); i++)
391  {
392  next_face_data.NodeIndices.push_back(mpVtkGeometryFilter->GetOutput()->GetCell(mBoundaryFacesRead + mBoundaryFacesSkipped)->GetPointId(i));
393  }
394  }
395  else if (ELEMENT_DIM == 2u)
396  {
397  for (unsigned i = 0; i < (mNodesPerElement-1); i++)
398  {
399  next_face_data.NodeIndices.push_back(mpVtkFilterEdges->GetOutput()->GetCell(mBoundaryFacesRead)->GetPointId(i));
400  }
401  }
402  else if (ELEMENT_DIM == 1u)
403  {
404  vtkSmartPointer<vtkIdList> enclosing_cells = vtkSmartPointer<vtkIdList>::New();
405  assert(mNumNodes == (unsigned) mpVtkUnstructuredGrid->GetNumberOfPoints());
406  for (unsigned point_index = mBoundaryFacesSkipped; point_index < mNumNodes; ++point_index)
407  {
408  mBoundaryFacesSkipped++;
409  mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
410 
411  if (enclosing_cells->GetNumberOfIds() == 1u)
412  {
413  next_face_data.NodeIndices.push_back(point_index);
414  break;
415  }
416  }
417  }
418  else
419  {
421  }
422 
423  mBoundaryFacesRead++;
424  return next_face_data;
425 }
426 
427 
428 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
429 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<double>& dataPayload)
430 {
431  vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
432 
433  if ( !p_cell_data->HasArray(dataName.c_str()) )
434  {
435  EXCEPTION("No cell data '" + dataName + "'");
436  }
437 
438  vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
439  if (p_scalars->GetNumberOfComponents() != 1)
440  {
441  EXCEPTION("The cell data '" + dataName + "' is not scalar data.");
442  }
443 
444  dataPayload.clear();
445  for (unsigned i = 0; i < mNumElements; i++)
446  {
447  dataPayload.push_back( p_scalars->GetTuple(i)[0] );
448  }
449 }
450 
451 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
452 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
453 {
454  vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
455 
456  if ( !p_cell_data->HasArray(dataName.c_str()) )
457  {
458  EXCEPTION("No cell data '" + dataName + "'");
459  }
460 
461  vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
462  if (p_scalars->GetNumberOfComponents() != 3)
463  {
464  EXCEPTION("The cell data '" + dataName + "' is not 3-vector data.");
465  }
466 
467  dataPayload.clear();
468  for (unsigned i = 0; i < mNumElements; i++)
469  {
470  c_vector <double, SPACE_DIM> data;
471  for (unsigned j = 0; j < SPACE_DIM; j++)
472  {
473  data[j]= p_scalars->GetTuple(i)[j];
474  }
475  dataPayload.push_back(data);
476  }
477 }
478 
479 
480 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
481 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<double>& dataPayload)
482 {
483  vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
484 
485  if ( !p_point_data->HasArray(dataName.c_str()) )
486  {
487  EXCEPTION("No point data '" + dataName + "'");
488  }
489 
490  vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
491  if (p_scalars->GetNumberOfComponents() != 1)
492  {
493  EXCEPTION("The point data '" + dataName + "' is not scalar data.");
494  }
495 
496  dataPayload.clear();
497 
498  for (unsigned i = 0; i < mNumNodes; i++)
499  {
500  dataPayload.push_back( p_scalars->GetTuple(i)[0] );
501  }
502 }
503 
504 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
505 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
506 {
507  vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
508 
509  if ( !p_point_data->HasArray(dataName.c_str()) )
510  {
511  EXCEPTION("No point data '" + dataName + "'");
512  }
513 
514  vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
515 
516  if (p_scalars->GetNumberOfComponents() != 3)
517  {
518  EXCEPTION("The point data '" + dataName + "' is not 3-vector data.");
519  }
520  dataPayload.clear();
521  for (unsigned i = 0; i < mNumNodes; i++)
522  {
523  c_vector<double, SPACE_DIM> data;
524  for (unsigned j = 0; j< SPACE_DIM; j++)
525  {
526  data[j] = p_scalars->GetTuple(i)[j];
527  }
528  dataPayload.push_back( data );
529  }
530 }
531 
532 
533 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
535 {
536  return mpVtkUnstructuredGrid;
537 }
538 
540 // Explicit instantiation
542 
543 template class VtkMeshReader<1,1>;
544 template class VtkMeshReader<1,2>;
545 template class VtkMeshReader<1,3>;
546 template class VtkMeshReader<2,2>;
547 template class VtkMeshReader<2,3>;
548 template class VtkMeshReader<3,3>;
549 #endif // CHASTE_VTK
void CommonConstructor()
unsigned GetNumCableElementAttributes() const
std::ifstream mVtuFile
#define EXCEPTION(message)
Definition: Exception.hpp:143
unsigned GetNumEdges() const
void GetPointData(std::string dataName, std::vector< double > &dataPayload)
unsigned GetNumFaceAttributes() const
ElementData GetNextFaceData()
unsigned GetNumElementAttributes() const
#define NEVER_REACHED
Definition: Exception.hpp:206
ElementData GetNextElementData()
std::vector< unsigned > NodeIndices
unsigned GetNumElements() const
virtual ~VtkMeshReader()
ElementData GetNextCableElementData()
unsigned GetNumCableElements() const
void GetCellData(std::string dataName, std::vector< double > &dataPayload)
unsigned GetNumFaces() const
unsigned GetNumNodes() const
vtkUnstructuredGrid * OutputMeshAsVtkUnstructuredGrid()
std::vector< double > GetNextNode()
vtkSmartPointer< vtkUnstructuredGrid > mpVtkUnstructuredGrid
VtkMeshReader(std::string pathBaseName)