]> Creatis software - clitk.git/blob - tools/clitkMeshToBinaryImage.cxx
Merge branch 'master' of git.creatis.insa-lyon.fr:clitk
[clitk.git] / tools / clitkMeshToBinaryImage.cxx
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to:
5   - University of LYON              http://www.universite-lyon.fr/
6   - Léon Bérard cancer center       http://www.centreleonberard.fr
7   - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
8
9   This software is distributed WITHOUT ANY WARRANTY; without even
10   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11   PURPOSE.  See the copyright notices for more information.
12
13   It is distributed under dual licence
14
15   - BSD        See included LICENSE.txt file
16   - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17 ===========================================================================*/
18 #include "clitkMeshToBinaryImage_ggo.h"
19 #include "clitkMeshToBinaryImageFilter.h"
20 #include "clitkCommon.h"
21 #include "clitkImageCommon.h"
22
23 #include <itkImageIOBase.h>
24 #include <itkInvertIntensityImageFilter.h>
25
26 #include <vtkOBJReader.h>
27
28 #include <string>
29
30 template <unsigned dim>
31 void run(const args_info_clitkMeshToBinaryImage& argsInfo,
32          const itk::ImageIOBase * header)
33 {
34   typedef itk::Image<unsigned char, dim> ImageType;
35   typename ImageType::Pointer like = ImageType::New();
36
37   typename ImageType::SizeType size;
38   for (unsigned i = 0; i < dim; ++i)
39     size[i] = header->GetDimensions(i);
40   typename ImageType::RegionType region;
41   region.SetSize(size);
42   like->SetRegions(region);
43
44   typename ImageType::SpacingType spacing;
45   for (unsigned i = 0; i < dim; ++i)
46     spacing[i] = header->GetSpacing(i);
47   like->SetSpacing(spacing);
48
49
50   typename ImageType::PointType origin;
51   for (unsigned i = 0; i < dim; ++i)
52     origin[i] = header->GetOrigin(i);
53   like->SetOrigin(origin);
54
55   vtkSmartPointer<vtkOBJReader> reader = vtkOBJReader::New();
56   reader->SetFileName(argsInfo.input_arg);
57   reader->Update();
58
59   typename clitk::MeshToBinaryImageFilter<ImageType>::Pointer filter =
60     clitk::MeshToBinaryImageFilter<ImageType>::New();
61   filter->SetExtrude(false);
62   filter->SetMesh(reader->GetOutput());
63   filter->SetLikeImage(like);
64   filter->Update();
65
66   typedef itk::InvertIntensityImageFilter<ImageType> InvertFilterType;
67   typename InvertFilterType::Pointer ifilter = InvertFilterType::New();
68   ifilter->SetInput(filter->GetOutput());
69   ifilter->SetMaximum(1);
70
71   clitk::writeImage(ifilter->GetOutput(), argsInfo.output_arg);
72 }
73
74 int main(int argc, char** argv)
75 {
76   GGO(clitkMeshToBinaryImage, args_info);
77
78   itk::ImageIOBase::Pointer header =
79     clitk::readImageHeader(args_info.like_arg);
80   switch (header->GetNumberOfDimensions())
81   {
82     case 2:
83       run<2>(args_info, header);
84       break;
85     case 3:
86       run<3>(args_info, header);
87       break;
88     case 4:
89       run<4>(args_info, header);
90       break;
91   }
92
93   return EXIT_SUCCESS;
94 }