Chaste Release::3.1
VtkMeshReader.cpp
00001 /*
00002 
00003 Copyright (C) Fujitsu Laboratories of Europe, 2009
00004 
00005 */
00006 
00007 /*
00008 
00009 Copyright (c) 2005-2012, University of Oxford.
00010 All rights reserved.
00011 
00012 University of Oxford means the Chancellor, Masters and Scholars of the
00013 University of Oxford, having an administrative office at Wellington
00014 Square, Oxford OX1 2JD, UK.
00015 
00016 This file is part of Chaste.
00017 
00018 Redistribution and use in source and binary forms, with or without
00019 modification, are permitted provided that the following conditions are met:
00020  * Redistributions of source code must retain the above copyright notice,
00021    this list of conditions and the following disclaimer.
00022  * Redistributions in binary form must reproduce the above copyright notice,
00023    this list of conditions and the following disclaimer in the documentation
00024    and/or other materials provided with the distribution.
00025  * Neither the name of the University of Oxford nor the names of its
00026    contributors may be used to endorse or promote products derived from this
00027    software without specific prior written permission.
00028 
00029 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00030 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00031 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00032 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00033 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00034 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00035 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00036 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00037 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00038 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
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     // Check file exists
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     // Load the mesh geometry and data from a file
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 == 2)
00096     {
00097         mNodesPerElement = 3;
00098         mVtkCellType = VTK_TRIANGLE;
00099     }
00100 
00101     //Determine if we have multiple cell types - such as cable elements in addition to tets/triangles
00102     vtkCellTypes* cell_types = vtkCellTypes::New();
00103     mpVtkUnstructuredGrid->GetCellTypes(cell_types);
00104 
00105     if(cell_types->GetNumberOfTypes() > 1)
00106     {
00107         mNumCableElementAttributes = 1;
00108         for(unsigned cell_id = 0; cell_id < num_cells; ++cell_id)
00109         {
00110             if(mpVtkUnstructuredGrid->GetCellType(cell_id) == mVtkCellType)
00111             {
00112                 ++mNumElements;
00113                 assert(mNumCableElements == 0); //We expect all the simplices first and then the cables at the end of the array
00114             }
00115             else if(mpVtkUnstructuredGrid->GetCellType(cell_id) == VTK_LINE)
00116             {
00117                 ++mNumCableElements;
00118             }
00119             else
00120             {
00121                 NEVER_REACHED;
00122             }
00123         }
00124     }
00125     else
00126     {
00127         //There is only 1 cell type, so all cells are elements
00128         mNumElements = num_cells;
00129     }
00130 
00131     cell_types->Delete();
00132 
00133 
00134     // Extract the surface faces
00135     if (SPACE_DIM == 2)
00136     {
00137         vtkDataSetSurfaceFilter* p_surface = vtkDataSetSurfaceFilter::New();
00138         p_surface->SetInput(mpVtkUnstructuredGrid);
00139         mpVtkFilterEdges = vtkFeatureEdges::New();
00140         mpVtkFilterEdges->SetInput(p_surface->GetOutput());
00141         mpVtkFilterEdges->Update();
00142         mNumFaces = mpVtkFilterEdges->GetOutput()->GetNumberOfCells();
00143         p_surface->Delete();
00144     }
00145     else
00146     {
00147         mpVtkGeometryFilter = vtkGeometryFilter::New();
00148         mpVtkGeometryFilter->SetInput(mpVtkUnstructuredGrid);
00149         mpVtkGeometryFilter->Update();
00150 
00151         mNumFaces = mpVtkGeometryFilter->GetOutput()->GetNumberOfCells();
00152         if (mNumCableElements > 0)
00153         {
00154             //The boundary face filter includes the cable elements - get rid of them
00155             unsigned num_cells = mNumFaces;
00156             for (unsigned i=0; i<num_cells; i++)
00157             {
00158                if (mpVtkGeometryFilter->GetOutput()->GetCellType(i) == VTK_LINE)
00159                {
00160                    mNumFaces--;
00161                }
00162             }
00163         }
00164 
00165     }
00166 }
00167 
00168 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00169 VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::VtkMeshReader(vtkUnstructuredGrid* p_vtkUnstructuredGrid) :
00170     mIndexFromZero(true),
00171     mNumNodes(0),
00172     mNumElements(0),
00173     mNumFaces(0),
00174     mNumCableElements(0),
00175     mNodesRead(0),
00176     mElementsRead(0),
00177     mFacesRead(0),
00178     mBoundaryFacesRead(0),
00179     mBoundaryFacesSkipped(0),
00180     mCableElementsRead(0),
00181     mNumElementAttributes(0),
00182     mNumFaceAttributes(0),
00183     mNumCableElementAttributes(0),
00184     mOrderOfElements(1),
00185     mNodesPerElement(4),
00186     mVtkCellType(VTK_TETRA)
00187 {
00188     mpVtkUnstructuredGrid = p_vtkUnstructuredGrid;
00189     CommonConstructor();
00190 }
00191 
00192 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00193 VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::~VtkMeshReader()
00194 {
00195     if (SPACE_DIM == 3)
00196     {
00197         mpVtkGeometryFilter->Delete();
00198     }
00199     else
00200     {
00201         mpVtkFilterEdges->Delete();
00202     }
00203 }
00204 
00205 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00206 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumElements() const
00207 {
00208     return mNumElements;
00209 }
00210 
00211 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00212 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumCableElements() const
00213 {
00214     return mNumCableElements;
00215 }
00216 
00217 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00218 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumNodes() const
00219 {
00220     return mNumNodes;
00221 }
00222 
00223 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00224 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumFaces() const
00225 {
00226     return mNumFaces;
00227 }
00228 
00229 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00230 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumEdges() const
00231 {
00232     return mNumFaces;
00233 }
00234 
00235 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00236 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumElementAttributes() const
00237 {
00238     return mNumElementAttributes;
00239 }
00240 
00241 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00242 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumCableElementAttributes() const
00243 {
00244     return mNumCableElementAttributes;
00245 }
00246 
00247 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00248 unsigned VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNumFaceAttributes() const
00249 {
00250     return mNumFaceAttributes;
00251 }
00252 
00253 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00254 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::Reset()
00255 {
00256     mNodesRead=0;
00257     mElementsRead=0;
00258     mFacesRead=0;
00259     mBoundaryFacesRead=0;
00260     mBoundaryFacesSkipped=0;
00261     mCableElementsRead=0;
00262 }
00263 
00264 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00265 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::Initialize()
00266 {
00267     mpVtkUnstructuredGrid->Initialize();
00268 }
00269 
00270 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00271 std::vector<double> VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNextNode()
00272 {
00273     if ( mNodesRead >= mNumNodes )
00274     {
00275         EXCEPTION( "Trying to read data for a node that doesn't exist" );
00276     }
00277 
00278     std::vector<double> next_node;
00279 
00280     for (unsigned i = 0; i < 3; i++)
00281     {
00282         next_node.push_back( mpVtkUnstructuredGrid->GetPoint(mNodesRead)[i] );
00283     }
00284 
00285     mNodesRead++;
00286     return next_node;
00287 }
00288 
00289 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00290 ElementData VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNextElementData()
00291 {
00292     if ( mElementsRead >= mNumElements )
00293     {
00294         EXCEPTION( "Trying to read data for an element that doesn't exist" );
00295     }
00296 
00297     if (mpVtkUnstructuredGrid->GetCellType(mElementsRead)!=  mVtkCellType)
00298     {
00299         EXCEPTION("Element is not of expected type (vtkTetra/vtkTriangle)");
00300     }
00301 
00302     ElementData next_element_data;
00303 
00304     for (unsigned i = 0; i < mNodesPerElement; i++)
00305     {
00306         next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(mElementsRead)->GetPointId(i));
00307     }
00308 
00309     // \todo implement method to read element data properly (currently returns zero always...)
00310     next_element_data.AttributeValue = 0;
00311 
00312     mElementsRead++;
00313     return next_element_data;
00314 }
00315 
00316 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00317 ElementData VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNextCableElementData()
00318 {
00319     if ( mCableElementsRead >=  mNumCableElements )
00320     {
00321         EXCEPTION( "Trying to read data for a cable element that doesn't exist" );
00322     }
00323 
00324     unsigned next_index = mNumElements + mCableElementsRead;
00325     assert(mpVtkUnstructuredGrid->GetCellType(next_index)==VTK_LINE);
00326 
00327     ElementData next_element_data;
00328 
00329     for (unsigned i = 0; i < 2; i++)
00330     {
00331         next_element_data.NodeIndices.push_back(mpVtkUnstructuredGrid->GetCell(next_index)->GetPointId(i));
00332     }
00333 
00334     vtkDataArray *p_scalars = mpVtkUnstructuredGrid->GetCellData()->GetArray( "Cable radius" );
00335     next_element_data.AttributeValue = p_scalars->GetTuple(next_index)[0];
00336 
00337     mCableElementsRead++;
00338     return next_element_data;
00339 }
00340 
00341 
00342 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00343 ElementData VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetNextFaceData()
00344 {
00345     if ( mBoundaryFacesRead >= mNumFaces)
00346     {
00347         EXCEPTION( "Trying to read data for a boundary element that doesn't exist");
00348     }
00349 
00350     ElementData next_face_data;
00351 
00352     if (SPACE_DIM == 3)
00353     {
00354         while (mpVtkGeometryFilter->GetOutput()->GetCellType(mBoundaryFacesRead + mBoundaryFacesSkipped) == VTK_LINE)
00355         {
00356             mBoundaryFacesSkipped++;
00357         }
00358         for (unsigned i = 0; i < (mNodesPerElement-1); i++)
00359         {
00360             next_face_data.NodeIndices.push_back(mpVtkGeometryFilter->GetOutput()->GetCell(mBoundaryFacesRead + mBoundaryFacesSkipped)->GetPointId(i));
00361         }
00362     }
00363     else
00364     {
00365         for (unsigned i = 0; i < (mNodesPerElement-1); i++)
00366         {
00367             next_face_data.NodeIndices.push_back(mpVtkFilterEdges->GetOutput()->GetCell(mBoundaryFacesRead)->GetPointId(i));
00368         }
00369     }
00370     mBoundaryFacesRead++;
00371     return next_face_data;
00372 }
00373 
00374 
00375 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00376 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<double>& dataPayload)
00377 {
00378     vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
00379 
00380     if ( !p_cell_data->HasArray(dataName.c_str()) )
00381     {
00382         EXCEPTION("No cell data '" + dataName + "'");
00383     }
00384 
00385     vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
00386     if (p_scalars->GetNumberOfComponents() != 1)
00387     {
00388         EXCEPTION("The cell data '" + dataName + "' is not scalar data.");
00389     }
00390 
00391     dataPayload.clear();
00392     for (unsigned i = 0; i < mNumElements; i++)
00393     {
00394         dataPayload.push_back( p_scalars->GetTuple(i)[0] );
00395     }
00396 }
00397 
00398 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00399 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetCellData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
00400 {
00401     vtkCellData *p_cell_data = mpVtkUnstructuredGrid->GetCellData();
00402 
00403     if ( !p_cell_data->HasArray(dataName.c_str()) )
00404     {
00405         EXCEPTION("No cell data '" + dataName + "'");
00406     }
00407 
00408     vtkDataArray *p_scalars = p_cell_data->GetArray( dataName.c_str() );
00409     if (p_scalars->GetNumberOfComponents() != 3)
00410     {
00411         EXCEPTION("The cell data '" + dataName + "' is not 3-vector data.");
00412     }
00413 
00414     dataPayload.clear();
00415     for (unsigned i = 0; i < mNumElements; i++)
00416     {
00417         c_vector <double, SPACE_DIM> data;
00418         for (unsigned j = 0; j < SPACE_DIM; j++)
00419         {
00420             data[j]=  p_scalars->GetTuple(i)[j];
00421         }
00422         dataPayload.push_back(data);
00423     }
00424 }
00425 
00426 
00427 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00428 void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::GetPointData(std::string dataName, std::vector<double>& dataPayload)
00429 {
00430     vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
00431 
00432     if ( !p_point_data->HasArray(dataName.c_str()) )
00433     {
00434         EXCEPTION("No point data '" + dataName + "'");
00435     }
00436 
00437     vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
00438     if (p_scalars->GetNumberOfComponents() != 1)
00439     {
00440         EXCEPTION("The point data '" + dataName + "' is not scalar data.");
00441     }
00442 
00443     dataPayload.clear();
00444 
00445     for (unsigned i = 0; i < mNumNodes; 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>::GetPointData(std::string dataName, std::vector<c_vector<double,SPACE_DIM> >& dataPayload)
00453 {
00454    vtkPointData *p_point_data = mpVtkUnstructuredGrid->GetPointData();
00455 
00456     if ( !p_point_data->HasArray(dataName.c_str()) )
00457     {
00458         EXCEPTION("No point data '" + dataName + "'");
00459     }
00460 
00461     vtkDataArray *p_scalars = p_point_data->GetArray( dataName.c_str() );
00462 
00463     if (p_scalars->GetNumberOfComponents() != 3)
00464     {
00465         EXCEPTION("The point data '" + dataName + "' is not 3-vector data.");
00466     }
00467     dataPayload.clear();
00468     for (unsigned i = 0; i < mNumNodes; 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 vtkUnstructuredGrid* VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::OutputMeshAsVtkUnstructuredGrid()
00482 {
00483     return mpVtkUnstructuredGrid;
00484 }
00485 
00487 // Explicit instantiation
00489 
00490 template class VtkMeshReader<1,1>;
00491 template class VtkMeshReader<1,2>;
00492 template class VtkMeshReader<1,3>;
00493 template class VtkMeshReader<2,2>;
00494 template class VtkMeshReader<2,3>;
00495 template class VtkMeshReader<3,3>;
00496 #endif // CHASTE_VTK