1 #include <tclap/CmdLine.h> 3 #include <itkMaskImageFilter.h> 4 #include <itkVectorImage.h> 6 #include <itkImageRegionIterator.h> 7 #include <itkExtractImageFilter.h> 18 template <
class ImageType>
22 using MaskImageType = itk::Image <unsigned short, 3>;
23 MaskImageType::Pointer maskData = anima::readImage <MaskImageType> (args.
maskName);
25 if (args.
imageIO->GetNumberOfComponents() > 1)
26 throw itk::ExceptionObject(__FILE__, __LINE__,
"4D vector images are not supported yet", ITK_LOCATION);
28 using InternalImageType = itk::Image <typename ImageType::PixelType, 3>;
29 using MaskImageFilterType = itk::MaskImageFilter <InternalImageType, MaskImageType, InternalImageType>;
31 typename ImageType::Pointer dataImage = anima::readImage <ImageType> (args.
refName);
32 unsigned int size4d = dataImage->GetLargestPossibleRegion().GetSize()[3];
34 using ExtractFilterType = itk::ExtractImageFilter <ImageType, InternalImageType>;
35 for (
unsigned int i = 0;i < size4d;++i)
37 typename ImageType::RegionType region = dataImage->GetLargestPossibleRegion();
41 typename ExtractFilterType::Pointer extractFilter = ExtractFilterType::New();
42 extractFilter->SetInput(dataImage);
43 extractFilter->SetExtractionRegion(region);
44 extractFilter->SetDirectionCollapseToGuess();
46 extractFilter->Update();
48 typename MaskImageFilterType::Pointer maskFilter = MaskImageFilterType::New();
50 maskFilter->SetInput(extractFilter->GetOutput());
51 maskFilter->SetMaskImage(maskData);
54 using ImageIteratorType = itk::ImageRegionIterator <ImageType>;
55 using InternalIteratorType = itk::ImageRegionIterator <InternalImageType>;
58 ImageIteratorType dataItr(dataImage, region);
60 InternalIteratorType maskedItr(maskFilter->GetOutput(), extractFilter->GetOutput()->GetLargestPossibleRegion());
61 while (!maskedItr.IsAtEnd())
63 dataItr.Set(maskedItr.Get());
69 anima::writeImage <ImageType> (args.
resName, dataImage);
72 template <
class ImageType>
76 using MaskImageType = itk::Image <unsigned short, 3>;
77 MaskImageType::Pointer maskData = anima::readImage <MaskImageType> (args.
maskName);
79 using MaskImageFilterType = itk::MaskImageFilter <ImageType, MaskImageType, ImageType>;
80 typename MaskImageFilterType::Pointer maskFilter = MaskImageFilterType::New();
82 maskFilter->SetInput(anima::readImage <ImageType> (args.
refName));
83 maskFilter->SetMaskImage(maskData);
86 anima::writeImage <ImageType> (args.
resName, maskFilter->GetOutput());
89 template <
class ComponentType,
int dimension>
93 if (imageIO->GetNumberOfComponents() > 1)
96 throw itk::ExceptionObject (__FILE__, __LINE__,
"Number of dimensions not supported for vector image masking", ITK_LOCATION);
98 mask3DImage < itk::VectorImage<ComponentType, 3> > (args);
103 mask3DImage < itk::Image<ComponentType, 3> > (args);
105 maskScalar4DImage < itk::Image<ComponentType, 4> > (args);
109 template <
class ComponentType>
116 int main(
int argc,
char **argv)
118 TCLAP::CmdLine cmd(
"INRIA / IRISA - VisAGeS/Empenn Team",
' ',ANIMA_VERSION);
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);
126 cmd.parse(argc,argv);
128 catch (TCLAP::ArgException& e)
130 std::cerr <<
"Error: " << e.error() <<
"for argument " << e.argId() << std::endl;
135 itk::ImageIOBase::Pointer
imageIO = itk::ImageIOFactory::CreateImageIO(inArg.getValue().c_str(),
136 itk::ImageIOFactory::ReadMode);
139 std::cerr <<
"Itk could not find suitable IO factory for the input " << inArg.getValue() << std::endl;
144 imageIO->SetFileName(inArg.getValue());
145 imageIO->ReadImageInformation();
148 args.
refName = inArg.getValue();
150 args.
resName = outArg.getValue();
157 catch (itk::ExceptionObject &err)
159 std::cerr <<
"Cannot perform maskin, be sure to use valid arguments..." << std::endl;
160 std::cerr << err << std::endl;
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,...)