Chaste Commit::675f9facbe008c5eacb9006feaeb6423206579ea
VtkMeshReader.cpp
1/*
2
3Copyright (C) Fujitsu Laboratories of Europe, 2009
4
5*/
6
7/*
8
9Copyright (c) 2005-2025, University of Oxford.
10All rights reserved.
11
12University of Oxford means the Chancellor, Masters and Scholars of the
13University of Oxford, having an administrative office at Wellington
14Square, Oxford OX1 2JD, UK.
15
16This file is part of Chaste.
17
18Redistribution and use in source and binary forms, with or without
19modification, 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
29THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
30AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
31IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
32ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
33LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
34CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
35GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
36HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
37LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
38OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
39
40*/
41
42#ifdef CHASTE_VTK
43
44#include <vtkCellTypes.h>
45
46#include "VtkMeshReader.hpp"
47#include "Exception.hpp"
48
49template <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
88template <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#if (VTK_MAJOR_VERSION >= 9 && VTK_MINOR_VERSION >= 2)
107 const auto num_distinct_cell_types = static_cast<unsigned>(mpVtkUnstructuredGrid->GetDistinctCellTypesArray()->GetNumberOfTuples());
108#else // VTK older than 9.2
109 vtkCellTypes* cell_types = vtkCellTypes::New();
110 mpVtkUnstructuredGrid->GetCellTypes(cell_types);
111 const auto num_distinct_cell_types = static_cast<unsigned>(cell_types->GetNumberOfTypes());
112 cell_types->Delete();
113#endif
114
115 if (num_distinct_cell_types > 1u)
116 {
117 mNumCableElementAttributes = 1;
118 for (unsigned cell_id = 0; cell_id < num_cells; ++cell_id)
119 {
120 if (mpVtkUnstructuredGrid->GetCellType(cell_id) == mVtkCellType)
121 {
122 ++mNumElements;
123 assert(mNumCableElements == 0); //We expect all the simplices first and then the cables at the end of the array
124 }
125 else if (mpVtkUnstructuredGrid->GetCellType(cell_id) == VTK_LINE)
126 {
127 ++mNumCableElements;
128 }
129 else
130 {
132 }
133 }
134 }
135 else
136 {
137 //There is only 1 cell type, so all cells are elements
138 mNumElements = num_cells;
139 }
140
141 // Extract the surface faces
142 if (ELEMENT_DIM == 2u)
143 {
144 vtkDataSetSurfaceFilter* p_surface = vtkDataSetSurfaceFilter::New();
145 mpVtkFilterEdges = vtkFeatureEdges::New();
146#if VTK_MAJOR_VERSION >= 6
147 p_surface->SetInputData(mpVtkUnstructuredGrid);
148 mpVtkFilterEdges->SetInputConnection(p_surface->GetOutputPort());
149#else
150 p_surface->SetInput(mpVtkUnstructuredGrid);
151 mpVtkFilterEdges->SetInput(p_surface->GetOutput());
152#endif
153 mpVtkFilterEdges->Update();
154 mNumFaces = mpVtkFilterEdges->GetOutput()->GetNumberOfCells();
155 p_surface->Delete();
156 }
157 else if (ELEMENT_DIM == 3u)
158 {
159 mpVtkGeometryFilter = vtkGeometryFilter::New();
160#if VTK_MAJOR_VERSION >= 6
161 mpVtkGeometryFilter->SetInputData(mpVtkUnstructuredGrid);
162#else
163 mpVtkGeometryFilter->SetInput(mpVtkUnstructuredGrid);
164#endif
165
166#if (VTK_MAJOR_VERSION >= 9 && VTK_MINOR_VERSION >= 1)
167 // Change to indexing in vtkGeometryFilter happened in VTK 9.1
168 mpVtkGeometryFilter->SetPassThroughPointIds(1);
169#endif
170 mpVtkGeometryFilter->Update();
171
172 mNumFaces = mpVtkGeometryFilter->GetOutput()->GetNumberOfCells();
173 if (mNumCableElements > 0)
174 {
175 //The boundary face filter includes the cable elements - get rid of them
176 unsigned num_all_cells = mNumFaces;
177 for (unsigned i=0; i<num_all_cells; i++)
178 {
179 if (mpVtkGeometryFilter->GetOutput()->GetCellType(i) == VTK_LINE)
180 {
181 mNumFaces--;
182 }
183 }
184 }
185 }
186 else if (ELEMENT_DIM == 1u)
187 {
188 mNumFaces = 0;
189 vtkSmartPointer<vtkIdList> enclosing_cells = vtkSmartPointer<vtkIdList>::New();
190 assert(mNumNodes == (unsigned) mpVtkUnstructuredGrid->GetNumberOfPoints());
191 for (unsigned point_index = 0; point_index < mNumNodes; ++point_index)
192 {
193 mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
194
195 if (enclosing_cells->GetNumberOfIds() == 1u)
196 {
197 mNumFaces++;
198 }
199 }
200 }
201 else
202 {
204 }
205}
206
207template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
208VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::VtkMeshReader(vtkUnstructuredGrid* p_vtkUnstructuredGrid) :
209 mIndexFromZero(true),
210 mNumNodes(0),
211 mNumElements(0),
212 mNumFaces(0),
213 mNumCableElements(0),
214 mNodesRead(0),
215 mElementsRead(0),
216 mFacesRead(0),
217 mBoundaryFacesRead(0),
218 mBoundaryFacesSkipped(0),
219 mCableElementsRead(0),
220 mNumElementAttributes(0),
221 mNumFaceAttributes(0),
222 mNumCableElementAttributes(0),
223 mOrderOfElements(1),
224 mNodesPerElement(4),
225 mVtkCellType(VTK_TETRA)
226{
227 mpVtkUnstructuredGrid = p_vtkUnstructuredGrid;
229}
230
231template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
233{
234 if (ELEMENT_DIM == 3u)
235 {
236 mpVtkGeometryFilter->Delete();
237 }
238 else if (ELEMENT_DIM == 2u)
239 {
240 mpVtkFilterEdges->Delete();
241 }
242}
243
244template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
246{
247 return mNumElements;
248}
249
250template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
252{
253 return mNumCableElements;
254}
255
256template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
258{
259 return mNumNodes;
260}
261
262template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
264{
265 return mNumFaces;
266}
267
268template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
270{
271 return mNumFaces;
272}
273
274template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
276{
277 return mNumElementAttributes;
278}
279
280template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
282{
283 return mNumCableElementAttributes;
284}
285
286template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
288{
289 return mNumFaceAttributes;
290}
291
292template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
294{
295 mNodesRead = 0;
296 mElementsRead = 0;
297 mFacesRead = 0;
298 mBoundaryFacesRead = 0;
299 mBoundaryFacesSkipped = 0;
300 mCableElementsRead = 0;
301}
302
303template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
305{
306 if (mNodesRead >= mNumNodes)
307 {
308 EXCEPTION( "Trying to read data for a node that doesn't exist" );
309 }
310
311 std::vector<double> next_node;
312
313 for (unsigned i = 0; i < 3; i++)
314 {
315 next_node.push_back( mpVtkUnstructuredGrid->GetPoint(mNodesRead)[i] );
316 }
317
318 mNodesRead++;
319 return next_node;
320}
321
322template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
324{
325 if (mElementsRead >= mNumElements)
326 {
327 EXCEPTION( "Trying to read data for an element that doesn't exist" );
328 }
329
330 if (mpVtkUnstructuredGrid->GetCellType(mElementsRead)!= mVtkCellType)
331 {
332 EXCEPTION("Element is not of expected type (vtkTetra/vtkTriangle)");
333 }
334
335 ElementData next_element_data;
336
337 for (unsigned i = 0; i < mNodesPerElement; i++)
338 {
339 next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(mElementsRead)->GetPointId(i));
340 }
341
342 // \todo implement method to read element data properly (currently returns zero always...)
343 next_element_data.AttributeValue = 0;
344
345 mElementsRead++;
346 return next_element_data;
347}
348
349template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
351{
352 if (mCableElementsRead >= mNumCableElements)
353 {
354 EXCEPTION( "Trying to read data for a cable element that doesn't exist" );
355 }
356
357 unsigned next_index = mNumElements + mCableElementsRead;
358 assert(mpVtkUnstructuredGrid->GetCellType(next_index)==VTK_LINE);
359
360 ElementData next_element_data;
361
362 for (unsigned i = 0; i < 2; i++)
363 {
364 next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(next_index)->GetPointId(i));
365 }
366
367 vtkDataArray *p_scalars = mpVtkUnstructuredGrid->GetCellData()->GetArray( "Cable radius" );
368 next_element_data.AttributeValue = p_scalars->GetTuple(next_index)[0];
369
370 mCableElementsRead++;
371 return next_element_data;
372}
373
374
375template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
377{
378 if (mBoundaryFacesRead >= mNumFaces)
379 {
380 EXCEPTION( "Trying to read data for a boundary element that doesn't exist");
381 }
382
383 ElementData next_face_data;
384 if (ELEMENT_DIM == 3u)
385 {
386#if (VTK_MAJOR_VERSION >= 9 && VTK_MINOR_VERSION >= 1)
387 // Change to indexing in vtkGeometryFilter happened in VTK 9.1
388 vtkDataArray *original_ids = mpVtkGeometryFilter->GetOutput()->GetPointData()->GetArray("vtkOriginalPointIds");
389#endif
390 while (mpVtkGeometryFilter->GetOutput()->GetCellType(mBoundaryFacesRead + mBoundaryFacesSkipped) == VTK_LINE)
391 {
392 mBoundaryFacesSkipped++;
393 }
394 for (unsigned i = 0; i < (mNodesPerElement-1); i++)
395 {
396 unsigned id = mpVtkGeometryFilter->GetOutput()->GetCell(mBoundaryFacesRead + mBoundaryFacesSkipped)->GetPointId(i);
397#if (VTK_MAJOR_VERSION >= 9 && VTK_MINOR_VERSION >= 1)
398 // Change to indexing in vtkGeometryFilter happened in VTK 9.1
399 next_face_data.NodeIndices.push_back(original_ids->GetTuple1(id));
400#else
401 next_face_data.NodeIndices.push_back(id);
402#endif
403 }
404 }
405 else if (ELEMENT_DIM == 2u)
406 {
407 for (unsigned i = 0; i < (mNodesPerElement-1); i++)
408 {
409 next_face_data.NodeIndices.push_back(mpVtkFilterEdges->GetOutput()->GetCell(mBoundaryFacesRead)->GetPointId(i));
410 }
411 }
412 else if (ELEMENT_DIM == 1u)
413 {
414 vtkSmartPointer<vtkIdList> enclosing_cells = vtkSmartPointer<vtkIdList>::New();
415 assert(mNumNodes == (unsigned) mpVtkUnstructuredGrid->GetNumberOfPoints());
416 for (unsigned point_index = mBoundaryFacesSkipped; point_index < mNumNodes; ++point_index)
417 {
418 mBoundaryFacesSkipped++;
419 mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
420
421 if (enclosing_cells->GetNumberOfIds() == 1u)
422 {
423 next_face_data.NodeIndices.push_back(point_index);
424 break;
425 }
426 }
427 }
428 else
429 {
431 }
432
433 mBoundaryFacesRead++;
434 return next_face_data;
435}
436
437
438template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
439void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<double>& dataPayload)
440{
441 vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
442
443 if (!p_cell_data->HasArray(dataName.c_str()))
444 {
445 EXCEPTION("No cell data '" + dataName + "'");
446 }
447
448 vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
449 if (p_scalars->GetNumberOfComponents() != 1)
450 {
451 EXCEPTION("The cell data '" + dataName + "' is not scalar data.");
452 }
453
454 dataPayload.clear();
455 for (unsigned i = 0; i < mNumElements; i++)
456 {
457 dataPayload.push_back( p_scalars->GetTuple(i)[0] );
458 }
459}
460
461template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
462void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
463{
464 vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
465
466 if (!p_cell_data->HasArray(dataName.c_str()))
467 {
468 EXCEPTION("No cell data '" + dataName + "'");
469 }
470
471 vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
472 if (p_scalars->GetNumberOfComponents() != 3)
473 {
474 EXCEPTION("The cell data '" + dataName + "' is not 3-vector data.");
475 }
476
477 dataPayload.clear();
478 for (unsigned i = 0; i < mNumElements; i++)
479 {
480 c_vector <double, SPACE_DIM> data;
481 for (unsigned j = 0; j < SPACE_DIM; j++)
482 {
483 data[j]= p_scalars->GetTuple(i)[j];
484 }
485 dataPayload.push_back(data);
486 }
487}
488
489
490template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
491void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<double>& dataPayload)
492{
493 vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
494
495 if (!p_point_data->HasArray(dataName.c_str()))
496 {
497 EXCEPTION("No point data '" + dataName + "'");
498 }
499
500 vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
501 if (p_scalars->GetNumberOfComponents() != 1)
502 {
503 EXCEPTION("The point data '" + dataName + "' is not scalar data.");
504 }
505
506 dataPayload.clear();
507
508 for (unsigned i = 0; i < mNumNodes; i++)
509 {
510 dataPayload.push_back( p_scalars->GetTuple(i)[0] );
511 }
512}
513
514template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
515void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
516{
517 vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
518
519 if (!p_point_data->HasArray(dataName.c_str()))
520 {
521 EXCEPTION("No point data '" + dataName + "'");
522 }
523
524 vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
525
526 if (p_scalars->GetNumberOfComponents() != 3)
527 {
528 EXCEPTION("The point data '" + dataName + "' is not 3-vector data.");
529 }
530 dataPayload.clear();
531 for (unsigned i = 0; i < mNumNodes; i++)
532 {
533 c_vector<double, SPACE_DIM> data;
534 for (unsigned j = 0; j< SPACE_DIM; j++)
535 {
536 data[j] = p_scalars->GetTuple(i)[j];
537 }
538 dataPayload.push_back( data );
539 }
540}
541
542
543template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
545{
546 return mpVtkUnstructuredGrid;
547}
548
549// Explicit instantiation
550template class VtkMeshReader<1,1>;
551template class VtkMeshReader<1,2>;
552template class VtkMeshReader<1,3>;
553template class VtkMeshReader<2,2>;
554template class VtkMeshReader<2,3>;
555template class VtkMeshReader<3,3>;
556#endif // CHASTE_VTK
#define EXCEPTION(message)
#define NEVER_REACHED
ElementData GetNextCableElementData()
unsigned GetNumCableElementAttributes() const
vtkSmartPointer< vtkUnstructuredGrid > mpVtkUnstructuredGrid
void GetPointData(std::string dataName, std::vector< double > &dataPayload)
unsigned GetNumElementAttributes() const
unsigned GetNumCableElements() const
void GetCellData(std::string dataName, std::vector< double > &dataPayload)
VtkMeshReader(std::string pathBaseName)
void CommonConstructor()
unsigned GetNumFaceAttributes() const
std::ifstream mVtuFile
std::vector< double > GetNextNode()
virtual ~VtkMeshReader()
vtkUnstructuredGrid * OutputMeshAsVtkUnstructuredGrid()
unsigned GetNumNodes() const
unsigned GetNumEdges() const
ElementData GetNextFaceData()
unsigned GetNumFaces() const
unsigned GetNumElements() const
ElementData GetNextElementData()
std::vector< unsigned > NodeIndices