[ITK] [ITK-users] Problems with mesh I/O and conversions in ITK

Carlos Henrique Villa Pinto chvillap at gmail.com
Mon Apr 14 14:49:40 EDT 2014


Hi Matt,

As far as I could see, every model from the brain database has information
about points, triangle strips and normals in their .vtk files. So using
visualization tools like Paraview and Slicer, I can see the (original)
models surfaces entirely, not just their points.

However, the outputs provided by the MeshFileWriter class keep only the
points' data, even though my program performs only a simple read-write.

At first I thought this problem could be due to some incompatibility of
file format versions, since the original models are VTK DataFiles v3.0 and
ITK saves the outputs as v2.0. But for that other model taken from the ITK
example, everything works fine despite this discrepancy.

Thanks,
[]s


2014-04-14 13:41 GMT-03:00 Matt McCormick <matt.mccormick at kitware.com>:

> Hi Carlos,
>
> When using MeshFileReader/MeshFileWriter, which cell data is being
> lost?  Examining a few of the files in the link provided at
> nac-brain-atlas-1.0/models/ I am only seeing points.
>
> Thanks,
> Matt
>
> On Mon, Apr 14, 2014 at 8:54 AM, Carlos Henrique Villa Pinto
> <chvillap at gmail.com> wrote:
> > Hi everyone!
> > I'm having some problems working with meshes in ITK.
> >
> > What I'm trying to do (by now) is just a brief test in which I read a
> .vtk
> > mesh file from this database:
> > http://www.spl.harvard.edu/publications/item/view/2037, write it again
> to
> > the disk and see if everything is ok. Very simple.
> >
> > The problem is: if I try to do that using the itk::MeshReader and
> > itk::MeshWriter classes, the program doesn't work as expected since the
> > written mesh apparently loses all information about the original mesh's
> > cells (I'm using Paraview and Slicer for visualization). The same thing
> > happens if I use the itk::VTKPolyDataReader and itk::VTKPolyDataWriter
> > classes for the mesh I/O. However, if I do the same test reading the .vtk
> > mesh file provided by this example:
> > http://itk.org/ITKExamples/src/IO/Mesh/ReadMesh/Documentation.html,
> > everything works fine in both cases. So I guess it could be a file format
> > issue.
> >
> > I also tried to use VTK instead of ITK for the mesh I/O (the
> > vtkPolyDataReader and vtkPolyDataWriter classes). And with that, the
> program
> > works fine even with the meshes from the brain database. But as I intend
> to
> > work with deformable models algorithms in ITK, that approach requires me
> to
> > convert the read vtkPolyData into a ITK simplex mesh, but such conversion
> > makes my program crash (either an exception is thrown because of an
> invalid
> > point or the program enters in an infinite loop). I'm using the
> > vtkPolyDataToitkMesh class (from InsightApplications) to convert the VTK
> > polydata to a ITK triangle mesh, and the
> > itk::TriangleMeshToSimplexMeshFilter class to convert the latter to a
> > simplex mesh.
> >
> > Last, but not least: if I just do the following conversions: VTK
> polydata ->
> > ITK triangle mesh -> VTK polydata and try to write the resulting polydata
> > (using vtkPolyDataWriter), I get a segmentation fault. No clue about what
> > causes it.
> >
> > The more recent and complete version of my code (VTK + ITK) is below. If
> > anyone could help me with any of these issues, I would be very grateful.
> >
> >
> =============================================================================================
> >
> >   // Set the ITK type definitions.
> >
> >   const int Dimension = 3;
> >   typedef double PixelType;
> >
> >   typedef vtkPolyDataToitkMesh::TriangleMeshType
> > TriangleMeshType;
> >   typedef itk::DefaultDynamicMeshTraits<PixelType, Dimension, Dimension,
> > double, double> SimplexMeshTraits;
> >   typedef itk::SimplexMesh<PixelType, Dimension, SimplexMeshTraits>
> > SimplexMeshType;
> >
> >   typedef itk::TriangleMeshToSimplexMeshFilter<TriangleMeshType,
> > SimplexMeshType> TriangleToSimplexFilterType;
> >   typedef itk::SimplexMeshToTriangleMeshFilter<SimplexMeshType,
> > TriangleMeshType> SimplexToTriangleFilterType;
> >
> >   typedef itk::SimplexMeshVolumeCalculator<SimplexMeshType>
> > VolumeCalculatorType;
> >
> >   typedef itk::VTKPolyDataWriter<TriangleMeshType> WriterType;
> >
> >   // Read the input VTK polydata.
> >
> >   vtkSmartPointer<vtkPolyDataReader> reader =
> > vtkSmartPointer<vtkPolyDataReader>::New();
> >   reader->SetFileName(inputFileName.c_str());
> >   reader->Update();
> >
> >   std::cout << "Polydata read ok\n";
> >
> >   // Convert the VTK polydata to a ITK mesh.
> >
> >   vtkPolyDataToitkMesh polyDataToMesh;
> >   polyDataToMesh.SetInput(reader->GetOutput());
> >   polyDataToMesh.ConvertvtkToitk();
> >
> >   std::cout << "Conversion to triangle mesh ok\n";
> >
> >   // Convert the triangle mesh to a simplex mesh.
> >
> >   TriangleToSimplexFilterType::Pointer triangleToSimplexFilter =
> > TriangleToSimplexFilterType::New();
> >   triangleToSimplexFilter->SetInput(polyDataToMesh.GetOutput());
> >   triangleToSimplexFilter->Update();
> >
> >   std::cout << "Conversion to simplex mesh ok\n";
> >
> >   // Compute the simplex mesh volume.
> >
> >   VolumeCalculatorType::Pointer volumeCalculator =
> > VolumeCalculatorType::New();
> >   volumeCalculator->SetSimplexMesh(triangleToSimplexFilter->GetOutput());
> >   volumeCalculator->Compute();
> >
> >   std::cout << "Volume calculation ok (value = " <<
> > volumeCalculator->GetVolume() << ")\n";
> >
> >   // Convert the simplex mesh to triangle mesh.
> >
> >   SimplexToTriangleFilterType::Pointer simplexToTriangleFilter =
> > SimplexToTriangleFilterType::New();
> >
> simplexToTriangleFilter->SetInput(triangleToSimplexFilter->GetOutput());
> >   simplexToTriangleFilter->Update();
> >
> >   std::cout << "Conversion to triangle mesh ok\n";
> >
> >   // Convert the ITK mesh to VTK polydata.
> >
> >   itkMeshTovtkPolyData meshToPolyData;
> >   meshToPolyData.SetInput(simplexToTriangleFilter->GetOutput());
> >   meshToPolyData.ConvertitkTovtk();
> >
> >   std::cout << "Conversion to polydata ok\n";
> >
> >   // Write the output polydata.
> >
> >   vtkSmartPointer<vtkPolyDataWriter> writer =
> > vtkSmartPointer<vtkPolyDataWriter>::New();
> >   writer->SetFileName(outputFileName.c_str());
> >   writer->SetInputData(meshToPolyData.GetOutput());
> >   writer->Update();
> >
> >   std::cout << "Polydata write ok\n";
> >
> >
> =============================================================================================
> >
> > Thanks in advance.
> > []s
> >
> > --
> > Carlos Henrique Villa Pinto
> > Graduate Student in Computer Science
> > Federal University of São Carlos - Brazil
> > XCS
> >
> > _____________________________________
> > Powered by www.kitware.com
> >
> > Visit other Kitware open-source projects at
> > http://www.kitware.com/opensource/opensource.html
> >
> > Kitware offers ITK Training Courses, for more information visit:
> > http://www.kitware.com/products/protraining.php
> >
> > Please keep messages on-topic and check the ITK FAQ at:
> > http://www.itk.org/Wiki/ITK_FAQ
> >
> > Follow this link to subscribe/unsubscribe:
> > http://www.itk.org/mailman/listinfo/insight-users
> >
>



-- 
Carlos Henrique Villa Pinto
Graduate Student in Computer Science
Federal University of São Carlos - Brazil
XCS
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://public.kitware.com/pipermail/community/attachments/20140414/9800839f/attachment-0002.html>
-------------- next part --------------
_____________________________________
Powered by www.kitware.com

Visit other Kitware open-source projects at
http://www.kitware.com/opensource/opensource.html

Kitware offers ITK Training Courses, for more information visit:
http://www.kitware.com/products/protraining.php

Please keep messages on-topic and check the ITK FAQ at:
http://www.itk.org/Wiki/ITK_FAQ

Follow this link to subscribe/unsubscribe:
http://www.itk.org/mailman/listinfo/insight-users


More information about the Community mailing list