[Insight-users] How to properly use itkTriangleMeshToImageFilter?
Amardeep Singh
amar.singh at gmx.de
Fri Jul 3 13:38:23 EDT 2009
Dear Luis
Thank you very much for your quick reply.
What I actually want to do is iterate over a part of my image and
extract a circular neighborhood
at every point and postprocess these extracted circles. Using the
itkMesh was a first try to get
the circular neighborhood. However, I think it might be too slow for my
task as it would be applied very frequently.
So, I need to come up with something else, anyway, I guess. Any ideas?
I added the following line to my code:
reader->GetOutput()->TransformIndexToPhysicalPoint(idx, center);
where "idx" is the Index in my image where I want my sphere to be and
"center" is the center of the sphere (RegularSphereMeshSource).
Like this, I manage to place the sphere where I want it to be. I was not
aware that the center is set in world coordinates.
Thanks a lot for your help!
Best regards
Amardeep
Luis Ibanez wrote:
>
> Hi Amardeep,
>
> I think that this typically happens when the image parameters
> are such that the image is not overlapping with your input mesh
> in Physical space.
>
>
> You should choose the
>
> * Origin
> * Spacing
> * Number of pixels
>
> of you image to create a grid on top of the Physical space that
> is occupied by the Mesh.
>
> Have you looked at your Mesh, and found the bounds of the
> mesh nodes in space ?
>
> Can you please post to the mailing list what the BoundingBox
> of your mesh is ?
>
> These values can be used to select proper parameters for the
> output image that contains the rasterization of the Mesh.
>
>
> The image parameters should match the extent of the input Mesh,
> not the extent of an image that you read from disk. Unless your
> Mesh was extracted from that image that you are reading from
> disk...
>
>
>
> You could get a feeling for these parameters by loading your
> Mesh into ParaView (or any other visualization tool).
>
>
> Please let us know,
>
>
> Thanks
>
>
> Luis
>
>
> --------------------------------------------------------------------------------
> On Thu, Jul 2, 2009 at 3:34 PM, Amardeep Singh <amar.singh at gmx.de
> <mailto:amar.singh at gmx.de>> wrote:
>
> Dear ITK Users
>
> I would like to use the "itkTriangleMeshToImageFilter". In a first
> experiment, I just want to create
> a sphere using the "itkRegularSphereMeshSource", binarize it and
> save the output. I take the size, origin and spacing
> of the resulting image from an image that I read from the harddisk
> (just for convenience).
> Unfortunately, I get the following error:
>
> Start Processing!
> [-90, 126, -72] <--- origin
> [182, 218, 182] <--- size of image
> [1, 1, 1] <--- spacing of image
> [0, 0, 0] <--- index
> Exception caught !
>
> itk::ExceptionObject (0x81a3b00)
> Location: "void itk::TriangleMeshToBinaryImageFilter<TInputMesh,
> TOutputImage>::GenerateData() [with TInputMesh = itk::Mesh<float,
> 3u, itk::DefaultStaticMeshTraits<float, 3u, 3u, float, float,
> float> >, TOutputImage = itk::Image<unsigned char, 3u>]"
> File:
> /workspace/InsightToolkit-3.14.0/Code/BasicFilters/itkTriangleMeshToBinaryImageFilter.txx
> Line: 224
> Description: itk::ERROR:
> TriangleMeshToBinaryImageFilter(0x8141ec8): No Image Indices Found.
>
> Above, I have included output, indicating the origin, size,
> spacing and index of the largest region of the image that I read
> from the harddisk. These are the values that I pass to the
> "itkTriangleMeshToImageFilter", so they should be fine, I think.
> There has been a similar posting on the mailing list
> ("http://www.nabble.com/converting-mesh-to-binary-image-td21170657.html"),
> but the solution to the problem was not made explicit,
> unfortunately. Please, find my source code attached.
> Does anyone know how to resolve this problem?
> Thank you very much!
>
> Best regards
> Amardeep
>
>
> #include "itkImageRegionIterator.h"
> #include "itkImageFileReader.h"
> #include "itkImageFileWriter.h"
> #include "itkTriangleMeshToBinaryImageFilter.h"
> #include "itkDefaultStaticMeshTraits.h"
> #include "itkMesh.h"
> #include "itkRegularSphereMeshSource.h"
>
> #include "iostream"
>
>
> int main(int argc, char *argv[] )
> {
> std::cout << "Start Processing!" << std::endl;
>
> typedef unsigned char BinaryPixelType;
> typedef float PixelType;
> typedef itk::DefaultStaticMeshTraits<float, 3, 3,float,float>
> TriangleMeshTraits;
> typedef itk::Mesh<float,3, TriangleMeshTraits> TriangleMeshType;
> typedef itk::RegularSphereMeshSource<TriangleMeshType>
> SphereMeshSourceType;
> typedef SphereMeshSourceType::PointType PointType;
>
> typedef itk::Image<BinaryPixelType, 3> BinaryImageType;
> typedef itk::Image<PixelType, 3> ImageType;
>
> typedef itk::TriangleMeshToBinaryImageFilter<TriangleMeshType,
> BinaryImageType> TriangleMeshToBinaryImageFilterType;
>
>
> TriangleMeshToBinaryImageFilterType::Pointer triangleMeshToImage
> = TriangleMeshToBinaryImageFilterType::New();
> SphereMeshSourceType::Pointer sphere = SphereMeshSourceType::New();
> sphere->SetScale(5);
> sphere->SetResolution(5);
>
> typedef itk::ImageFileReader<ImageType> ReaderType;
> typedef itk::ImageFileWriter<BinaryImageType> WriterType;
>
> ReaderType::Pointer reader = ReaderType::New();
> reader->SetFileName("/test_image.nii.gz");
> try
> {
> reader->Update();
> }
> catch( itk::ExceptionObject & excep )
> {
> std::cerr << "Exception caught !" << std::endl;
> std::cerr << excep << std::endl;
> }
>
> PointType center;
> center[0] = 10;
> center[1] = 10;
> center[2] = 10;
> sphere->SetCenter(center);
> triangleMeshToImage->SetInput(sphere->GetOutput());
>
> ImageType::Pointer orgImage = reader->GetOutput();
> const ImageType::PointType orgOrigin = orgImage->GetOrigin();
> const ImageType::SizeType orgSize =
> orgImage->GetLargestPossibleRegion().GetSize();
> const ImageType::SpacingType orgSpacing = orgImage->GetSpacing();
> const ImageType::IndexType orgIndex =
> orgImage->GetLargestPossibleRegion().GetIndex();
>
> std::cout << orgOrigin << std::endl;
> std::cout << orgSize << std::endl;
> std::cout << orgSpacing << std::endl;
> std::cout << orgIndex << std::endl;
>
> triangleMeshToImage->SetTolerance (1.0);
> triangleMeshToImage->SetSpacing (orgSpacing);
> triangleMeshToImage->SetOrigin(orgOrigin);
> triangleMeshToImage->SetSize (orgSize);
> triangleMeshToImage->SetIndex (orgIndex);
>
> WriterType::Pointer writer = WriterType::New();
>
> writer->SetFileName("/test_output.nii.gz");
> writer->SetInput(triangleMeshToImage->GetOutput());
>
> try
> {
> sphere->Update();
> triangleMeshToImage->Update();
> writer->Update();
> }
> catch( itk::ExceptionObject & excep )
> {
> std::cerr << "Exception caught !" << std::endl;
> std::cerr << excep << std::endl;
> }
>
> return 0;
>
> }
>
>
>
>
>
> _____________________________________
> Powered by www.kitware.com <http://www.kitware.com>
>
> Visit other Kitware open-source projects at
> http://www.kitware.com/opensource/opensource.html
>
> 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 Insight-users
mailing list