46 #include "VtkMeshReader.hpp"
49 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
59 mBoundaryFacesRead(0),
60 mBoundaryFacesSkipped(0),
61 mCableElementsRead(0),
62 mNumElementAttributes(0),
63 mNumFaceAttributes(0),
64 mNumCableElementAttributes(0),
67 mVtkCellType(VTK_TETRA)
69 vtkXMLUnstructuredGridReader* vtk_xml_unstructured_grid_reader;
75 EXCEPTION(
"Could not open VTU file: " + pathBaseName);
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();
85 vtk_xml_unstructured_grid_reader->Delete();
88 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
92 mNumNodes = mpVtkUnstructuredGrid->GetNumberOfPoints();
93 unsigned num_cells = mpVtkUnstructuredGrid->GetNumberOfCells();
95 if (ELEMENT_DIM == 2u)
98 mVtkCellType = VTK_TRIANGLE;
100 else if (ELEMENT_DIM == 1u)
102 mNodesPerElement = 2;
103 mVtkCellType = VTK_LINE;
107 vtkCellTypes* cell_types = vtkCellTypes::New();
108 mpVtkUnstructuredGrid->GetCellTypes(cell_types);
110 if(cell_types->GetNumberOfTypes() > 1)
112 mNumCableElementAttributes = 1;
113 for(
unsigned cell_id = 0; cell_id < num_cells; ++cell_id)
115 if(mpVtkUnstructuredGrid->GetCellType(cell_id) == mVtkCellType)
118 assert(mNumCableElements == 0);
120 else if(mpVtkUnstructuredGrid->GetCellType(cell_id) == VTK_LINE)
133 mNumElements = num_cells;
136 cell_types->Delete();
140 if (ELEMENT_DIM == 2u)
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());
148 p_surface->SetInput(mpVtkUnstructuredGrid);
149 mpVtkFilterEdges->SetInput(p_surface->GetOutput());
151 mpVtkFilterEdges->Update();
152 mNumFaces = mpVtkFilterEdges->GetOutput()->GetNumberOfCells();
155 else if (ELEMENT_DIM == 3u)
157 mpVtkGeometryFilter = vtkGeometryFilter::New();
158 #if VTK_MAJOR_VERSION >= 6
159 mpVtkGeometryFilter->SetInputData(mpVtkUnstructuredGrid);
161 mpVtkGeometryFilter->SetInput(mpVtkUnstructuredGrid);
163 mpVtkGeometryFilter->Update();
165 mNumFaces = mpVtkGeometryFilter->GetOutput()->GetNumberOfCells();
166 if (mNumCableElements > 0)
169 unsigned num_all_cells = mNumFaces;
170 for (
unsigned i=0; i<num_all_cells; i++)
172 if (mpVtkGeometryFilter->GetOutput()->GetCellType(i) == VTK_LINE)
179 else if (ELEMENT_DIM == 1u)
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)
186 mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
188 if (enclosing_cells->GetNumberOfIds() == 1u)
200 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
202 mIndexFromZero(true),
206 mNumCableElements(0),
210 mBoundaryFacesRead(0),
211 mBoundaryFacesSkipped(0),
212 mCableElementsRead(0),
213 mNumElementAttributes(0),
214 mNumFaceAttributes(0),
215 mNumCableElementAttributes(0),
218 mVtkCellType(VTK_TETRA)
224 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
227 if (ELEMENT_DIM == 3u)
229 mpVtkGeometryFilter->Delete();
231 else if (ELEMENT_DIM == 2u)
233 mpVtkFilterEdges->Delete();
237 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
243 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
246 return mNumCableElements;
249 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
255 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
261 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
267 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
270 return mNumElementAttributes;
273 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
276 return mNumCableElementAttributes;
279 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
282 return mNumFaceAttributes;
285 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
291 mBoundaryFacesRead=0;
292 mBoundaryFacesSkipped=0;
293 mCableElementsRead=0;
296 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
299 mpVtkUnstructuredGrid->Initialize();
302 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
305 if ( mNodesRead >= mNumNodes )
307 EXCEPTION(
"Trying to read data for a node that doesn't exist" );
310 std::vector<double> next_node;
312 for (
unsigned i = 0; i < 3; i++)
314 next_node.push_back( mpVtkUnstructuredGrid->GetPoint(mNodesRead)[i] );
321 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
324 if ( mElementsRead >= mNumElements )
326 EXCEPTION(
"Trying to read data for an element that doesn't exist" );
329 if (mpVtkUnstructuredGrid->GetCellType(mElementsRead)!= mVtkCellType)
331 EXCEPTION(
"Element is not of expected type (vtkTetra/vtkTriangle)");
336 for (
unsigned i = 0; i < mNodesPerElement; i++)
338 next_element_data.
NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(mElementsRead)->GetPointId(i));
345 return next_element_data;
348 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
351 if ( mCableElementsRead >= mNumCableElements )
353 EXCEPTION(
"Trying to read data for a cable element that doesn't exist" );
356 unsigned next_index = mNumElements + mCableElementsRead;
357 assert(mpVtkUnstructuredGrid->GetCellType(next_index)==VTK_LINE);
361 for (
unsigned i = 0; i < 2; i++)
363 next_element_data.
NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(next_index)->GetPointId(i));
366 vtkDataArray *p_scalars = mpVtkUnstructuredGrid->GetCellData()->GetArray(
"Cable radius" );
367 next_element_data.AttributeValue = p_scalars->GetTuple(next_index)[0];
369 mCableElementsRead++;
370 return next_element_data;
374 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
377 if ( mBoundaryFacesRead >= mNumFaces)
379 EXCEPTION(
"Trying to read data for a boundary element that doesn't exist");
384 if (ELEMENT_DIM == 3u)
386 while (mpVtkGeometryFilter->GetOutput()->GetCellType(mBoundaryFacesRead + mBoundaryFacesSkipped) == VTK_LINE)
388 mBoundaryFacesSkipped++;
390 for (
unsigned i = 0; i < (mNodesPerElement-1); i++)
392 next_face_data.
NodeIndices.push_back(mpVtkGeometryFilter->GetOutput()->GetCell(mBoundaryFacesRead + mBoundaryFacesSkipped)->GetPointId(i));
395 else if (ELEMENT_DIM == 2u)
397 for (
unsigned i = 0; i < (mNodesPerElement-1); i++)
399 next_face_data.
NodeIndices.push_back(mpVtkFilterEdges->GetOutput()->GetCell(mBoundaryFacesRead)->GetPointId(i));
402 else if (ELEMENT_DIM == 1u)
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)
408 mBoundaryFacesSkipped++;
409 mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
411 if (enclosing_cells->GetNumberOfIds() == 1u)
423 mBoundaryFacesRead++;
424 return next_face_data;
428 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
431 vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
433 if ( !p_cell_data->HasArray(dataName.c_str()) )
435 EXCEPTION(
"No cell data '" + dataName +
"'");
438 vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
439 if (p_scalars->GetNumberOfComponents() != 1)
441 EXCEPTION(
"The cell data '" + dataName +
"' is not scalar data.");
445 for (
unsigned i = 0; i < mNumElements; i++)
447 dataPayload.push_back( p_scalars->GetTuple(i)[0] );
451 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
454 vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
456 if ( !p_cell_data->HasArray(dataName.c_str()) )
458 EXCEPTION(
"No cell data '" + dataName +
"'");
461 vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
462 if (p_scalars->GetNumberOfComponents() != 3)
464 EXCEPTION(
"The cell data '" + dataName +
"' is not 3-vector data.");
468 for (
unsigned i = 0; i < mNumElements; i++)
470 c_vector <double, SPACE_DIM> data;
471 for (
unsigned j = 0; j < SPACE_DIM; j++)
473 data[j]= p_scalars->GetTuple(i)[j];
475 dataPayload.push_back(data);
480 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
483 vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
485 if ( !p_point_data->HasArray(dataName.c_str()) )
487 EXCEPTION(
"No point data '" + dataName +
"'");
490 vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
491 if (p_scalars->GetNumberOfComponents() != 1)
493 EXCEPTION(
"The point data '" + dataName +
"' is not scalar data.");
498 for (
unsigned i = 0; i < mNumNodes; i++)
500 dataPayload.push_back( p_scalars->GetTuple(i)[0] );
504 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
507 vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
509 if ( !p_point_data->HasArray(dataName.c_str()) )
511 EXCEPTION(
"No point data '" + dataName +
"'");
514 vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
516 if (p_scalars->GetNumberOfComponents() != 3)
518 EXCEPTION(
"The point data '" + dataName +
"' is not 3-vector data.");
521 for (
unsigned i = 0; i < mNumNodes; i++)
523 c_vector<double, SPACE_DIM> data;
524 for (
unsigned j = 0; j< SPACE_DIM; j++)
526 data[j] = p_scalars->GetTuple(i)[j];
528 dataPayload.push_back( data );
533 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
536 return mpVtkUnstructuredGrid;
unsigned GetNumCableElementAttributes() const
#define EXCEPTION(message)
unsigned GetNumEdges() const
void GetPointData(std::string dataName, std::vector< double > &dataPayload)
unsigned GetNumFaceAttributes() const
ElementData GetNextFaceData()
unsigned GetNumElementAttributes() const
ElementData GetNextElementData()
std::vector< unsigned > NodeIndices
unsigned GetNumElements() const
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)