]> Creatis software - clitk.git/commitdiff
With clitkImageCreate, copy the transformMatrix of the mhd file
authortbaudier <thomas.baudier@creatis.insa-lyon.fr>
Thu, 16 Mar 2017 14:48:54 +0000 (15:48 +0100)
committertbaudier <thomas.baudier@creatis.insa-lyon.fr>
Thu, 16 Mar 2017 14:48:54 +0000 (15:48 +0100)
tools/clitkImageCreate.cxx
tools/clitkImageCreate.ggo

index ab2bec514673ced6d5f1313d3477dfa94015feda..d84bbb51400c0d4a7fbfd45091baead16f4d4ce8 100644 (file)
 #include "clitkIO.h"
 
 template<class ImageType>
-void NewFilledImage(int * size, float * spacing, double * origin,
+void NewFilledImage(int * size, float * spacing, double * origin, double * direction,
                     double value,typename ImageType::Pointer output)
 {
   static const unsigned int Dim = ImageType::GetImageDimension();
+
   typename ImageType::SizeType mSize;
   mSize.Fill (0);
   for(unsigned int i=0; i<Dim; i++) mSize[i] = size[i];
+
   typename ImageType::RegionType mRegion;
   mRegion.SetSize(mSize);
+
   typename ImageType::SpacingType mSpacing;
   for(unsigned int i=0; i<Dim; i++) mSpacing[i] = spacing[i];
+
+  typename ImageType::DirectionType directionImage;
+  for(unsigned int i=0; i<Dim; i++)
+    for(unsigned int j=0; j<Dim; j++)
+      directionImage[i][j] = direction[i*Dim+j];
+
   output->SetRegions(mRegion);
   output->SetSpacing(mSpacing);
+  output->SetDirection(directionImage);
   output->Allocate();
+
   typename ImageType::PointType mOrigin;
   for(unsigned int i=0; i<Dim; i++) mOrigin[i] = origin[i];
   output->SetOrigin(mOrigin);
+
   typedef typename ImageType::PixelType PixelType;
   PixelType p = clitk::PixelTypeDownCast<double, PixelType>(value);
   output->FillBuffer(p);
@@ -75,11 +87,15 @@ int main(int argc, char * argv[])
     args_info.spacing_arg = new float[dim];
     args_info.origin_given = dim;
     args_info.origin_arg = new double[dim];
+    args_info.transformMatrix_given = dim*dim;
+    args_info.transformMatrix_arg = new double[dim*dim];
 
     for(int i=0; i<dim; i++) {
       args_info.size_arg[i] = header->GetDimensions(i);
       args_info.spacing_arg[i] = header->GetSpacing(i);
       args_info.origin_arg[i]=   header->GetOrigin(i);
+      for (int j=0; j<dim; ++j)
+        args_info.transformMatrix_arg[i*dim+j] = header->GetDirection(i)[j];
     }
   }
 
@@ -119,27 +135,50 @@ int main(int argc, char * argv[])
     for(int i=0; i<dim; i++) spacing[i] = args_info.spacing_arg[i];
   }
 
+  // direction
+  std::vector<double> direction;
+  direction.resize(dim*dim);
+  for(int i=0; i<dim; i++)
+  {
+    for(int j=0; j<dim; j++)
+    {
+      if (i == j)
+        direction[i*dim+j] = 1;
+      else
+        direction[i*dim+j] = 0;
+    }
+  }
+  if (args_info.transformMatrix_given) {
+    if (args_info.transformMatrix_given != dim*dim) {
+      std::cerr << "ERROR : please give the same number of values for --transfomMatrix and --spacing." << std::endl;
+      exit(-1);
+    }
+    for(int i=0; i<dim; i++)
+      for(int j=0; j<dim; j++)
+        direction[i*dim+j] = args_info.transformMatrix_arg[i*dim+j];
+  }
+
   // Create new image
   typedef float PixelType;
   if (dim == 2) {
     const int Dim=2;
     typedef itk::Image<PixelType, Dim> ImageType;
     ImageType::Pointer output = ImageType::New();
-    NewFilledImage<ImageType>(args_info.size_arg, &spacing[0], &origin[0], args_info.value_arg, output);
+    NewFilledImage<ImageType>(args_info.size_arg, &spacing[0], &origin[0], &direction[0], args_info.value_arg, output);
     clitk::writeImage<ImageType>(output, args_info.output_arg);
   }
   if (dim == 3) {
     const int Dim=3;
     typedef itk::Image<PixelType, Dim> ImageType;
     ImageType::Pointer output = ImageType::New();
-    NewFilledImage<ImageType>(args_info.size_arg, &spacing[0], &origin[0], args_info.value_arg, output);
+    NewFilledImage<ImageType>(args_info.size_arg, &spacing[0], &origin[0], &direction[0], args_info.value_arg, output);
     clitk::writeImage<ImageType>(output, args_info.output_arg);
   }
   if (dim == 4) {
     const int Dim=4;
     typedef itk::Image<PixelType, Dim> ImageType;
     ImageType::Pointer output = ImageType::New();
-    NewFilledImage<ImageType>(args_info.size_arg, &spacing[0], &origin[0], args_info.value_arg, output);
+    NewFilledImage<ImageType>(args_info.size_arg, &spacing[0], &origin[0], &direction[0], args_info.value_arg, output);
     clitk::writeImage<ImageType>(output, args_info.output_arg);
   }
 
index f197de74b44d3beb1262fa383fc5f458aff00d47..391801d346d925911ed78c75d764fb14a5ff159b 100644 (file)
@@ -11,6 +11,7 @@ option "like"         l "Size/spacing like this other image" string no
 option "size"          - "Number of pixels of each coordinate" int     no multiple
 option "spacing"       - "Spacing in mm between pixels"        float   no multiple
 option "origin"                - "Origin in mm"                        double  no multiple
+option "transformMatrix" - "Rotation matrix" double no multiple
 
 option "value"         - "Value for all voxels"  float  default="0.0"  no
 option "verbose"       v "Verbose"               flag   off