1 #include <tclap/CmdLine.h> 4 #include <itkVectorImage.h> 5 #include <itkCommand.h> 6 #include <itkExtractImageFilter.h> 17 template <
class ComponentType,
unsigned int InputDim>
22 typedef itk::VectorImage<ComponentType, InputDim> InputImageType;
23 typedef itk::Image<ComponentType, InputDim+1> OutputImageType;
25 typename InputImageType::Pointer inputImg = anima::readImage<InputImageType>(args.
input);
26 unsigned int nbComp = imageIO->GetNumberOfComponents();
28 std::cout <<
"number of dimension : " << InputDim << std::endl;
29 std::cout <<
"number of components : " << nbComp << std::endl;
32 typename OutputImageType::RegionType finalRegion;
33 typename OutputImageType::SizeType finalSize;
34 typename OutputImageType::IndexType finalIndex;
35 typename OutputImageType::PointType finalOrigin;
36 typename OutputImageType::SpacingType finalSpacing;
37 typename OutputImageType::DirectionType finalDirection;
39 for (
unsigned int d = 0; d < InputDim; ++d)
42 finalIndex[d] = inputImg->GetLargestPossibleRegion().GetIndex()[d];
43 finalSize[d] = inputImg->GetLargestPossibleRegion().GetSize()[d];
44 finalOrigin[d] = inputImg->GetOrigin()[d];
45 finalSpacing[d] = inputImg->GetSpacing()[d];
46 for(
unsigned int i = 0; i < InputDim; ++i)
47 finalDirection[d][i] = inputImg->GetDirection()[d][i];
49 finalIndex[InputDim] = 0;
50 finalSize[InputDim] = nbComp;
51 finalOrigin[InputDim] = args.
origin;
52 finalSpacing[InputDim] = args.
spacing;
54 for(
unsigned int i = 0; i < InputDim; ++i)
56 finalDirection[InputDim][i] = 0;
57 finalDirection[i][InputDim] = 0;
59 finalDirection[InputDim][InputDim] = 1;
61 finalRegion.SetIndex(finalIndex);
62 finalRegion.SetSize(finalSize);
64 typename OutputImageType::Pointer outputImg = OutputImageType::New();
65 outputImg->Initialize();
66 outputImg->SetRegions(finalRegion);
67 outputImg->SetOrigin(finalOrigin);
68 outputImg->SetSpacing(finalSpacing);
69 outputImg->SetDirection(finalDirection);
70 outputImg->Allocate();
72 typedef itk::ImageRegionIterator <OutputImageType> FillIteratorType;
73 FillIteratorType fillItr(outputImg, outputImg->GetLargestPossibleRegion());
74 typedef itk::ImageRegionConstIterator<InputImageType> SourceIteratorType;
75 SourceIteratorType srcItr(inputImg, inputImg->GetLargestPossibleRegion());
77 while(!fillItr.IsAtEnd())
85 fillItr.Set(srcItr.Get()[idx]);
89 anima::writeImage<OutputImageType>(args.
output, outputImg);
92 template <
class ComponentType >
100 int main(
int ac,
const char** av)
103 TCLAP::CmdLine cmd(
"Collapse a vector image into a scalar image of dim n+1.\n\ 104 The last dimension of the output has same size as the component of the input image.\n\ 105 You can give the spacing and the origin of the added dimension, default values are 0 for the origin and 1 for the spacing\n\ 106 Example: input image is 4x4x4 with a component size of 32, the output will be 4x4x4x32.\n\ 107 INRIA / IRISA - VisAGeS/Empenn Team",
110 TCLAP::ValueArg<std::string> inputArg(
"i",
112 "Input images to collapse",
118 TCLAP::ValueArg<std::string> outputArg(
"o",
125 TCLAP::ValueArg<double> origArg(
"O",
127 "Origin of the added Dimension image",
132 TCLAP::ValueArg<double> spacingArg(
"s",
134 "Spacing of the added Dimension image",
144 catch (TCLAP::ArgException& e)
146 std::cerr <<
"Error: " << e.error() <<
"for argument " << e.argId() << std::endl;
151 itk::ImageIOBase::Pointer
imageIO = itk::ImageIOFactory::CreateImageIO(inputArg.getValue().c_str(),
152 itk::ImageIOFactory::ReadMode);
156 std::cerr <<
"Itk could not find suitable IO factory for the input" << std::endl;
161 imageIO->SetFileName(inputArg.getValue());
162 imageIO->ReadImageInformation();
164 if(imageIO->GetNumberOfComponents() < 2)
166 std::cerr <<
"Input have to have component of size > 1. Nothing to do.";
170 std::cout<<
"\npreparing filter...\n";
174 args.
input = inputArg.getValue();
175 args.
output = outputArg.getValue();
176 args.
origin = origArg.getValue();
177 args.
spacing = spacingArg.getValue();
184 catch ( itk::ExceptionObject & err )
186 std::cerr <<
"Itk cannot collapse, be sure to use valid input..." << std::endl;
187 std::cerr << err << std::endl;
itk::ImageIOBase::Pointer imageIO
int main(int ac, const char **av)
void retrieveNbDimensions(itk::ImageIOBase::Pointer imageIO, const arguments &args)
#define ANIMA_RETRIEVE_COMPONENT_TYPE(imageIO, function,...)
void collapse(itk::ImageIOBase::Pointer imageIO, const arguments &args)
#define ANIMA_RETRIEVE_NUMBER_OF_DIMENSIONS(imageIO, ComponentType, function,...)