00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044 #ifdef CHASTE_VTK
00045
00046 #include "VtkMeshReader.hpp"
00047 #include "Exception.hpp"
00048
00049 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00050 VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::VtkMeshReader(std::string pathBaseName) :
00051 mIndexFromZero(true),
00052 mNumNodes(0),
00053 mNumElements(0),
00054 mNumFaces(0),
00055 mNumCableElements(0),
00056 mNodesRead(0),
00057 mElementsRead(0),
00058 mFacesRead(0),
00059 mBoundaryFacesRead(0),
00060 mBoundaryFacesSkipped(0),
00061 mCableElementsRead(0),
00062 mNumElementAttributes(0),
00063 mNumFaceAttributes(0),
00064 mNumCableElementAttributes(0),
00065 mOrderOfElements(1),
00066 mNodesPerElement(4),
00067 mVtkCellType(VTK_TETRA)
00068 {
00069 vtkXMLUnstructuredGridReader* vtk_xml_unstructured_grid_reader;
00070
00071
00072 mVtuFile.open(pathBaseName.c_str());
00073 if ( !mVtuFile.is_open() )
00074 {
00075 EXCEPTION("Could not open VTU file: " + pathBaseName);
00076 }
00077 mVtuFile.close();
00078
00079
00080 vtk_xml_unstructured_grid_reader = vtkXMLUnstructuredGridReader::New();
00081 vtk_xml_unstructured_grid_reader->SetFileName( pathBaseName.c_str() );
00082 vtk_xml_unstructured_grid_reader->Update();
00083 mpVtkUnstructuredGrid = vtk_xml_unstructured_grid_reader->GetOutput();
00084 CommonConstructor();
00085 vtk_xml_unstructured_grid_reader->Delete();
00086 }
00087
00088 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00089 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::CommonConstructor()
00090 {
00091
00092 mNumNodes = mpVtkUnstructuredGrid->GetNumberOfPoints();
00093 unsigned num_cells = mpVtkUnstructuredGrid->GetNumberOfCells();
00094
00095 if (ELEMENT_DIM == 2u)
00096 {
00097 mNodesPerElement = 3;
00098 mVtkCellType = VTK_TRIANGLE;
00099 }
00100 else if (ELEMENT_DIM == 1u)
00101 {
00102 mNodesPerElement = 2;
00103 mVtkCellType = VTK_LINE;
00104 }
00105
00106
00107 vtkCellTypes* cell_types = vtkCellTypes::New();
00108 mpVtkUnstructuredGrid->GetCellTypes(cell_types);
00109
00110 if(cell_types->GetNumberOfTypes() > 1)
00111 {
00112 mNumCableElementAttributes = 1;
00113 for(unsigned cell_id = 0; cell_id < num_cells; ++cell_id)
00114 {
00115 if(mpVtkUnstructuredGrid->GetCellType(cell_id) == mVtkCellType)
00116 {
00117 ++mNumElements;
00118 assert(mNumCableElements == 0);
00119 }
00120 else if(mpVtkUnstructuredGrid->GetCellType(cell_id) == VTK_LINE)
00121 {
00122 ++mNumCableElements;
00123 }
00124 else
00125 {
00126 NEVER_REACHED;
00127 }
00128 }
00129 }
00130 else
00131 {
00132
00133 mNumElements = num_cells;
00134 }
00135
00136 cell_types->Delete();
00137
00138
00139
00140 if (ELEMENT_DIM == 2u)
00141 {
00142 vtkDataSetSurfaceFilter* p_surface = vtkDataSetSurfaceFilter::New();
00143 mpVtkFilterEdges = vtkFeatureEdges::New();
00144 #if VTK_MAJOR_VERSION >= 6
00145 p_surface->SetInputData(mpVtkUnstructuredGrid);
00146 mpVtkFilterEdges->SetInputConnection(p_surface->GetOutputPort());
00147 #else
00148 p_surface->SetInput(mpVtkUnstructuredGrid);
00149 mpVtkFilterEdges->SetInput(p_surface->GetOutput());
00150 #endif
00151 mpVtkFilterEdges->Update();
00152 mNumFaces = mpVtkFilterEdges->GetOutput()->GetNumberOfCells();
00153 p_surface->Delete();
00154 }
00155 else if (ELEMENT_DIM == 3u)
00156 {
00157 mpVtkGeometryFilter = vtkGeometryFilter::New();
00158 #if VTK_MAJOR_VERSION >= 6
00159 mpVtkGeometryFilter->SetInputData(mpVtkUnstructuredGrid);
00160 #else
00161 mpVtkGeometryFilter->SetInput(mpVtkUnstructuredGrid);
00162 #endif
00163 mpVtkGeometryFilter->Update();
00164
00165 mNumFaces = mpVtkGeometryFilter->GetOutput()->GetNumberOfCells();
00166 if (mNumCableElements > 0)
00167 {
00168
00169 unsigned num_all_cells = mNumFaces;
00170 for (unsigned i=0; i<num_all_cells; i++)
00171 {
00172 if (mpVtkGeometryFilter->GetOutput()->GetCellType(i) == VTK_LINE)
00173 {
00174 mNumFaces--;
00175 }
00176 }
00177 }
00178 }
00179 else if (ELEMENT_DIM == 1u)
00180 {
00181 mNumFaces = 0;
00182 vtkSmartPointer<vtkIdList> enclosing_cells = vtkSmartPointer<vtkIdList>::New();
00183 assert(mNumNodes == (unsigned) mpVtkUnstructuredGrid->GetNumberOfPoints());
00184 for (unsigned point_index = 0; point_index < mNumNodes; ++point_index)
00185 {
00186 mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
00187
00188 if (enclosing_cells->GetNumberOfIds() == 1u)
00189 {
00190 mNumFaces++;
00191 }
00192 }
00193 }
00194 else
00195 {
00196 NEVER_REACHED;
00197 }
00198 }
00199
00200 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00201 VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::VtkMeshReader(vtkUnstructuredGrid* p_vtkUnstructuredGrid) :
00202 mIndexFromZero(true),
00203 mNumNodes(0),
00204 mNumElements(0),
00205 mNumFaces(0),
00206 mNumCableElements(0),
00207 mNodesRead(0),
00208 mElementsRead(0),
00209 mFacesRead(0),
00210 mBoundaryFacesRead(0),
00211 mBoundaryFacesSkipped(0),
00212 mCableElementsRead(0),
00213 mNumElementAttributes(0),
00214 mNumFaceAttributes(0),
00215 mNumCableElementAttributes(0),
00216 mOrderOfElements(1),
00217 mNodesPerElement(4),
00218 mVtkCellType(VTK_TETRA)
00219 {
00220 mpVtkUnstructuredGrid = p_vtkUnstructuredGrid;
00221 CommonConstructor();
00222 }
00223
00224 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00225 VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::~VtkMeshReader()
00226 {
00227 if (ELEMENT_DIM == 3u)
00228 {
00229 mpVtkGeometryFilter->Delete();
00230 }
00231 else if (ELEMENT_DIM == 2u)
00232 {
00233 mpVtkFilterEdges->Delete();
00234 }
00235 }
00236
00237 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00238 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumElements() const
00239 {
00240 return mNumElements;
00241 }
00242
00243 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00244 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumCableElements() const
00245 {
00246 return mNumCableElements;
00247 }
00248
00249 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00250 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumNodes() const
00251 {
00252 return mNumNodes;
00253 }
00254
00255 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00256 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumFaces() const
00257 {
00258 return mNumFaces;
00259 }
00260
00261 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00262 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumEdges() const
00263 {
00264 return mNumFaces;
00265 }
00266
00267 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00268 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumElementAttributes() const
00269 {
00270 return mNumElementAttributes;
00271 }
00272
00273 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00274 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumCableElementAttributes() const
00275 {
00276 return mNumCableElementAttributes;
00277 }
00278
00279 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00280 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumFaceAttributes() const
00281 {
00282 return mNumFaceAttributes;
00283 }
00284
00285 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00286 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::Reset()
00287 {
00288 mNodesRead=0;
00289 mElementsRead=0;
00290 mFacesRead=0;
00291 mBoundaryFacesRead=0;
00292 mBoundaryFacesSkipped=0;
00293 mCableElementsRead=0;
00294 }
00295
00296 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00297 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::Initialize()
00298 {
00299 mpVtkUnstructuredGrid->Initialize();
00300 }
00301
00302 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00303 std::vector<double> VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNextNode()
00304 {
00305 if ( mNodesRead >= mNumNodes )
00306 {
00307 EXCEPTION( "Trying to read data for a node that doesn't exist" );
00308 }
00309
00310 std::vector<double> next_node;
00311
00312 for (unsigned i = 0; i < 3; i++)
00313 {
00314 next_node.push_back( mpVtkUnstructuredGrid->GetPoint(mNodesRead)[i] );
00315 }
00316
00317 mNodesRead++;
00318 return next_node;
00319 }
00320
00321 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00322 ElementData VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNextElementData()
00323 {
00324 if ( mElementsRead >= mNumElements )
00325 {
00326 EXCEPTION( "Trying to read data for an element that doesn't exist" );
00327 }
00328
00329 if (mpVtkUnstructuredGrid->GetCellType(mElementsRead)!= mVtkCellType)
00330 {
00331 EXCEPTION("Element is not of expected type (vtkTetra/vtkTriangle)");
00332 }
00333
00334 ElementData next_element_data;
00335
00336 for (unsigned i = 0; i < mNodesPerElement; i++)
00337 {
00338 next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(mElementsRead)->GetPointId(i));
00339 }
00340
00341
00342 next_element_data.AttributeValue = 0;
00343
00344 mElementsRead++;
00345 return next_element_data;
00346 }
00347
00348 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00349 ElementData VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNextCableElementData()
00350 {
00351 if ( mCableElementsRead >= mNumCableElements )
00352 {
00353 EXCEPTION( "Trying to read data for a cable element that doesn't exist" );
00354 }
00355
00356 unsigned next_index = mNumElements + mCableElementsRead;
00357 assert(mpVtkUnstructuredGrid->GetCellType(next_index)==VTK_LINE);
00358
00359 ElementData next_element_data;
00360
00361 for (unsigned i = 0; i < 2; i++)
00362 {
00363 next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(next_index)->GetPointId(i));
00364 }
00365
00366 vtkDataArray *p_scalars = mpVtkUnstructuredGrid->GetCellData()->GetArray( "Cable radius" );
00367 next_element_data.AttributeValue = p_scalars->GetTuple(next_index)[0];
00368
00369 mCableElementsRead++;
00370 return next_element_data;
00371 }
00372
00373
00374 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00375 ElementData VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNextFaceData()
00376 {
00377 if ( mBoundaryFacesRead >= mNumFaces)
00378 {
00379 EXCEPTION( "Trying to read data for a boundary element that doesn't exist");
00380 }
00381
00382 ElementData next_face_data;
00383
00384 if (ELEMENT_DIM == 3u)
00385 {
00386 while (mpVtkGeometryFilter->GetOutput()->GetCellType(mBoundaryFacesRead + mBoundaryFacesSkipped) == VTK_LINE)
00387 {
00388 mBoundaryFacesSkipped++;
00389 }
00390 for (unsigned i = 0; i < (mNodesPerElement-1); i++)
00391 {
00392 next_face_data.NodeIndices.push_back(mpVtkGeometryFilter->GetOutput()->GetCell(mBoundaryFacesRead + mBoundaryFacesSkipped)->GetPointId(i));
00393 }
00394 }
00395 else if (ELEMENT_DIM == 2u)
00396 {
00397 for (unsigned i = 0; i < (mNodesPerElement-1); i++)
00398 {
00399 next_face_data.NodeIndices.push_back(mpVtkFilterEdges->GetOutput()->GetCell(mBoundaryFacesRead)->GetPointId(i));
00400 }
00401 }
00402 else if (ELEMENT_DIM == 1u)
00403 {
00404 vtkSmartPointer<vtkIdList> enclosing_cells = vtkSmartPointer<vtkIdList>::New();
00405 assert(mNumNodes == (unsigned) mpVtkUnstructuredGrid->GetNumberOfPoints());
00406 for (unsigned point_index = mBoundaryFacesSkipped; point_index < mNumNodes; ++point_index)
00407 {
00408 mBoundaryFacesSkipped++;
00409 mpVtkUnstructuredGrid->GetPointCells(point_index, enclosing_cells);
00410
00411 if (enclosing_cells->GetNumberOfIds() == 1u)
00412 {
00413 next_face_data.NodeIndices.push_back(point_index);
00414 break;
00415 }
00416 }
00417 }
00418 else
00419 {
00420 NEVER_REACHED;
00421 }
00422
00423 mBoundaryFacesRead++;
00424 return next_face_data;
00425 }
00426
00427
00428 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00429 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<double>& dataPayload)
00430 {
00431 vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
00432
00433 if ( !p_cell_data->HasArray(dataName.c_str()) )
00434 {
00435 EXCEPTION("No cell data '" + dataName + "'");
00436 }
00437
00438 vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
00439 if (p_scalars->GetNumberOfComponents() != 1)
00440 {
00441 EXCEPTION("The cell data '" + dataName + "' is not scalar data.");
00442 }
00443
00444 dataPayload.clear();
00445 for (unsigned i = 0; i < mNumElements; i++)
00446 {
00447 dataPayload.push_back( p_scalars->GetTuple(i)[0] );
00448 }
00449 }
00450
00451 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00452 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
00453 {
00454 vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
00455
00456 if ( !p_cell_data->HasArray(dataName.c_str()) )
00457 {
00458 EXCEPTION("No cell data '" + dataName + "'");
00459 }
00460
00461 vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
00462 if (p_scalars->GetNumberOfComponents() != 3)
00463 {
00464 EXCEPTION("The cell data '" + dataName + "' is not 3-vector data.");
00465 }
00466
00467 dataPayload.clear();
00468 for (unsigned i = 0; i < mNumElements; i++)
00469 {
00470 c_vector <double, SPACE_DIM> data;
00471 for (unsigned j = 0; j < SPACE_DIM; j++)
00472 {
00473 data[j]= p_scalars->GetTuple(i)[j];
00474 }
00475 dataPayload.push_back(data);
00476 }
00477 }
00478
00479
00480 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00481 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<double>& dataPayload)
00482 {
00483 vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
00484
00485 if ( !p_point_data->HasArray(dataName.c_str()) )
00486 {
00487 EXCEPTION("No point data '" + dataName + "'");
00488 }
00489
00490 vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
00491 if (p_scalars->GetNumberOfComponents() != 1)
00492 {
00493 EXCEPTION("The point data '" + dataName + "' is not scalar data.");
00494 }
00495
00496 dataPayload.clear();
00497
00498 for (unsigned i = 0; i < mNumNodes; i++)
00499 {
00500 dataPayload.push_back( p_scalars->GetTuple(i)[0] );
00501 }
00502 }
00503
00504 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00505 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
00506 {
00507 vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
00508
00509 if ( !p_point_data->HasArray(dataName.c_str()) )
00510 {
00511 EXCEPTION("No point data '" + dataName + "'");
00512 }
00513
00514 vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
00515
00516 if (p_scalars->GetNumberOfComponents() != 3)
00517 {
00518 EXCEPTION("The point data '" + dataName + "' is not 3-vector data.");
00519 }
00520 dataPayload.clear();
00521 for (unsigned i = 0; i < mNumNodes; i++)
00522 {
00523 c_vector<double, SPACE_DIM> data;
00524 for (unsigned j = 0; j< SPACE_DIM; j++)
00525 {
00526 data[j] = p_scalars->GetTuple(i)[j];
00527 }
00528 dataPayload.push_back( data );
00529 }
00530 }
00531
00532
00533 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00534 vtkUnstructuredGrid* VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::OutputMeshAsVtkUnstructuredGrid()
00535 {
00536 return mpVtkUnstructuredGrid;
00537 }
00538
00540
00542
00543 template class VtkMeshReader<1,1>;
00544 template class VtkMeshReader<1,2>;
00545 template class VtkMeshReader<1,3>;
00546 template class VtkMeshReader<2,2>;
00547 template class VtkMeshReader<2,3>;
00548 template class VtkMeshReader<3,3>;
00549 #endif // CHASTE_VTK