ANIMA  4.0
animaMaskImage.cxx
Go to the documentation of this file.
1 #include <tclap/CmdLine.h>
2 
3 #include <itkMaskImageFilter.h>
4 #include <itkVectorImage.h>
5 #include <itkImage.h>
6 #include <itkImageRegionIterator.h>
7 #include <itkExtractImageFilter.h>
8 
11 
12 struct arguments
13 {
14  std::string refName, maskName, resName;
15  itk::ImageIOBase::Pointer imageIO;
16 };
17 
18 template <class ImageType>
19 void
21 {
22  using MaskImageType = itk::Image <unsigned short, 3>;
23  MaskImageType::Pointer maskData = anima::readImage <MaskImageType> (args.maskName);
24 
25  if (args.imageIO->GetNumberOfComponents() > 1)
26  throw itk::ExceptionObject(__FILE__, __LINE__, "4D vector images are not supported yet", ITK_LOCATION);
27 
28  using InternalImageType = itk::Image <typename ImageType::PixelType, 3>;
29  using MaskImageFilterType = itk::MaskImageFilter <InternalImageType, MaskImageType, InternalImageType>;
30 
31  typename ImageType::Pointer dataImage = anima::readImage <ImageType> (args.refName);
32  unsigned int size4d = dataImage->GetLargestPossibleRegion().GetSize()[3];
33 
34  using ExtractFilterType = itk::ExtractImageFilter <ImageType, InternalImageType>;
35  for (unsigned int i = 0;i < size4d;++i)
36  {
37  typename ImageType::RegionType region = dataImage->GetLargestPossibleRegion();
38  region.SetIndex(3,i);
39  region.SetSize(3,0);
40 
41  typename ExtractFilterType::Pointer extractFilter = ExtractFilterType::New();
42  extractFilter->SetInput(dataImage);
43  extractFilter->SetExtractionRegion(region);
44  extractFilter->SetDirectionCollapseToGuess();
45 
46  extractFilter->Update();
47 
48  typename MaskImageFilterType::Pointer maskFilter = MaskImageFilterType::New();
49 
50  maskFilter->SetInput(extractFilter->GetOutput());
51  maskFilter->SetMaskImage(maskData);
52  maskFilter->Update();
53 
54  using ImageIteratorType = itk::ImageRegionIterator <ImageType>;
55  using InternalIteratorType = itk::ImageRegionIterator <InternalImageType>;
56 
57  region.SetSize(3,1);
58  ImageIteratorType dataItr(dataImage, region);
59 
60  InternalIteratorType maskedItr(maskFilter->GetOutput(), extractFilter->GetOutput()->GetLargestPossibleRegion());
61  while (!maskedItr.IsAtEnd())
62  {
63  dataItr.Set(maskedItr.Get());
64  ++dataItr;
65  ++maskedItr;
66  }
67  }
68 
69  anima::writeImage <ImageType> (args.resName, dataImage);
70 }
71 
72 template <class ImageType>
73 void
74 mask3DImage(const arguments &args)
75 {
76  using MaskImageType = itk::Image <unsigned short, 3>;
77  MaskImageType::Pointer maskData = anima::readImage <MaskImageType> (args.maskName);
78 
79  using MaskImageFilterType = itk::MaskImageFilter <ImageType, MaskImageType, ImageType>;
80  typename MaskImageFilterType::Pointer maskFilter = MaskImageFilterType::New();
81 
82  maskFilter->SetInput(anima::readImage <ImageType> (args.refName));
83  maskFilter->SetMaskImage(maskData);
84  maskFilter->Update();
85 
86  anima::writeImage <ImageType> (args.resName, maskFilter->GetOutput());
87 }
88 
89 template <class ComponentType, int dimension>
90 void
91 checkIfComponentsAreVectors(itk::ImageIOBase::Pointer imageIO, const arguments &args)
92 {
93  if (imageIO->GetNumberOfComponents() > 1)
94  {
95  if (dimension > 3)
96  throw itk::ExceptionObject (__FILE__, __LINE__, "Number of dimensions not supported for vector image masking", ITK_LOCATION);
97 
98  mask3DImage < itk::VectorImage<ComponentType, 3> > (args);
99  }
100  else
101  {
102  if (dimension < 4)
103  mask3DImage < itk::Image<ComponentType, 3> > (args);
104  else
105  maskScalar4DImage < itk::Image<ComponentType, 4> > (args);
106  }
107 }
108 
109 template <class ComponentType>
110 void
111 retrieveNbDimensions(itk::ImageIOBase::Pointer imageIO, const arguments &args)
112 {
113  ANIMA_RETRIEVE_NUMBER_OF_DIMENSIONS(imageIO, ComponentType, checkIfComponentsAreVectors, imageIO, args);
114 }
115 
116 int main(int argc, char **argv)
117 {
118  TCLAP::CmdLine cmd("INRIA / IRISA - VisAGeS/Empenn Team", ' ',ANIMA_VERSION);
119 
120  TCLAP::ValueArg<std::string> inArg("i","inputfile","Input image",true,"","input image",cmd);
121  TCLAP::ValueArg<std::string> outArg("o","outputfile","output image",true,"","output image",cmd);
122  TCLAP::ValueArg<std::string> maskArg("m","maskfile","mask file",true,"","mask file",cmd);
123 
124  try
125  {
126  cmd.parse(argc,argv);
127  }
128  catch (TCLAP::ArgException& e)
129  {
130  std::cerr << "Error: " << e.error() << "for argument " << e.argId() << std::endl;
131  return EXIT_FAILURE;
132  }
133 
134  // Find out the type of the image in file
135  itk::ImageIOBase::Pointer imageIO = itk::ImageIOFactory::CreateImageIO(inArg.getValue().c_str(),
136  itk::ImageIOFactory::ReadMode);
137  if (!imageIO)
138  {
139  std::cerr << "Itk could not find suitable IO factory for the input " << inArg.getValue() << std::endl;
140  return EXIT_FAILURE;
141  }
142 
143  // Now that we found the appropriate ImageIO class, ask it to read the meta data from the image file.
144  imageIO->SetFileName(inArg.getValue());
145  imageIO->ReadImageInformation();
146 
147  arguments args;
148  args.refName = inArg.getValue();
149  args.maskName = maskArg.getValue();
150  args.resName = outArg.getValue();
151  args.imageIO = imageIO;
152 
153  try
154  {
155  ANIMA_RETRIEVE_COMPONENT_TYPE(imageIO, retrieveNbDimensions, imageIO, args)
156  }
157  catch (itk::ExceptionObject &err)
158  {
159  std::cerr << "Cannot perform maskin, be sure to use valid arguments..." << std::endl;
160  std::cerr << err << std::endl;
161  return EXIT_FAILURE;
162  }
163 
164  return EXIT_SUCCESS;
165 }
std::string maskName
void retrieveNbDimensions(itk::ImageIOBase::Pointer imageIO, const arguments &args)
int main(int argc, char **argv)
void checkIfComponentsAreVectors(itk::ImageIOBase::Pointer imageIO, const arguments &args)
void maskScalar4DImage(const arguments &args)
itk::ImageIOBase::Pointer imageIO
void mask3DImage(const arguments &args)
#define ANIMA_RETRIEVE_COMPONENT_TYPE(imageIO, function,...)
#define ANIMA_RETRIEVE_NUMBER_OF_DIMENSIONS(imageIO, ComponentType, function,...)
std::string refName
std::string resName