DMRITool  v0.1.1-139-g860d86b4
Diffusion MRI Tool
mexITKImageWrite.cxx
Go to the documentation of this file.
1 
19 #include "mex.h"
20 #include "utlMEX.h"
21 
22 template <class ImageType>
23 void
24 SetITKImageInformation ( itk::SmartPointer<ImageType>& image, const mxArray* originArray, const mxArray* spacingArray )
25 {
26  typename ImageType::PointType origin;
27  for ( int i = 0; i < ImageType::ImageDimension; i += 1 )
28  origin[i]=0;
29  if (originArray)
30  {
31  mwSize dimOrgin = mxGetNumberOfElements(originArray);
32  double * data = mxGetPr(originArray);
33  for ( int i = 0; i < utl::min((int)dimOrgin, (int)ImageType::ImageDimension); i += 1 )
34  origin[i] = data[i];
35  }
36  typename ImageType::SpacingType spacing;
37  for ( int i = 0; i < ImageType::ImageDimension; i += 1 )
38  spacing[i]=1;
39  if (spacingArray)
40  {
41  mwSize dimSpacing = mxGetNumberOfElements(spacingArray);
42  double * data = mxGetPr(spacingArray);
43  for ( int i = 0; i < utl::min((int)dimSpacing, (int)ImageType::ImageDimension); i += 1 )
44  spacing[i] = data[i];
45  }
46  image->SetOrigin(origin);
47  image->SetSpacing(spacing);
48 }
49 
50 template <typename T>
51  inline void callFunction(mxArray* plhs[], const mxArray* prhs[],
52  const int nlhs,const int nrhs, const bool isVectorImage)
53 {
54  std::string filename = utl::GetString(prhs[0]);
55 
56  mxArray* originArray=NULL;
57  mxArray* spacingArray=NULL;
58  if (nrhs==3)
59  {
60  originArray = utl::GetArrayStruct(prhs[2], "origin" );
61  spacingArray = utl::GetArrayStruct(prhs[2], "spacing" );
62  }
63 
64  mwSize dimArray = mxGetNumberOfDimensions(prhs[1]);
65  const mwSize* dims = mxGetDimensions(prhs[1]);
66 
67  if (isVectorImage)
68  {
69  utlGlobalException(dimArray!=4, "The image has to have 4 dimension");
70  typedef itk::VectorImage<T, 3> ImageType;
71  typename ImageType::Pointer image = ImageType::New();
72  itk::GetITKVectorImageFromMXArray(prhs[1], image);
73  SetITKImageInformation(image, originArray, spacingArray);
74  itk::SaveImage(image, filename);
75  }
76  else
77  {
78  typedef itk::Image<T, 4> ImageType;
79  typename ImageType::Pointer image = ImageType::New();
80  itk::GetITKImageFromMXArray(prhs[1], image);
81  SetITKImageInformation(image, originArray, spacingArray);
82  itk::SaveImage(image, filename);
83  }
84 
85 }
86 
87 
88 void mexFunction(int nlhs, mxArray *plhs[],
89  int nrhs, const mxArray *prhs[])
90 {
91  utlGlobalException(nrhs!=2 && nrhs!=3, "Bad number of inputs arguments");
92  utlGlobalException(nlhs!=0, "the program should no output");
93 
94  utlGlobalException(mxGetClassID(prhs[0])!=mxCHAR_CLASS, "the first input has to be a string");
95 
96  bool isVectorImage = utl::GetScalarStructDef<bool>(nrhs==3?prhs[2]:NULL,"vectorImage",false);
97 
98  if (mxGetClassID(prhs[1]) == mxSINGLE_CLASS)
99  {
100  std::cout << "mxSINGLE_CLASS, float" << std::endl << std::flush;
101  callFunction<float>(plhs,prhs,nlhs,nrhs, isVectorImage);
102  }
103  else if (mxGetClassID(prhs[1]) == mxLOGICAL_CLASS)
104  {
105  std::cout << "mxLOGICAL_CLASS, mxLogical" << std::endl << std::flush;
106  callFunction<mxLogical>(plhs,prhs,nlhs,nrhs, isVectorImage);
107  }
108  else if (mxGetClassID(prhs[1]) == mxUINT8_CLASS)
109  {
110  std::cout << "mxUINT8_CLASS, unsigned char" << std::endl << std::flush;
111  callFunction<unsigned char>(plhs,prhs,nlhs,nrhs, isVectorImage);
112  }
113  else if (mxGetClassID(prhs[1]) == mxUINT16_CLASS)
114  {
115  std::cout << "mxUINT16_CLASS, unsigned short int" << std::endl << std::flush;
116  callFunction<unsigned short int>(plhs,prhs,nlhs,nrhs, isVectorImage);
117  }
118  else if (mxGetClassID(prhs[1]) == mxUINT32_CLASS)
119  {
120  std::cout << "mxUINT32_CLASS, unsigned int" << std::endl << std::flush;
121  callFunction<unsigned int>(plhs,prhs,nlhs,nrhs, isVectorImage);
122  }
123  else if (mxGetClassID(prhs[1]) == mxUINT64_CLASS)
124  {
125  std::cout << "mxUINT64_CLASS, uint64_T" << std::endl << std::flush;
126  callFunction<uint64_T>(plhs,prhs,nlhs,nrhs, isVectorImage);
127  }
128  else if (mxGetClassID(prhs[1]) == mxINT8_CLASS)
129  {
130  std::cout << "mxINT8_CLASS, signed char" << std::endl << std::flush;
131  callFunction<signed char>(plhs,prhs,nlhs,nrhs, isVectorImage);
132  }
133  else if (mxGetClassID(prhs[1]) == mxINT16_CLASS)
134  {
135  std::cout << "mxINT16_CLASS, short int" << std::endl << std::flush;
136  callFunction<short int>(plhs,prhs,nlhs,nrhs, isVectorImage);
137  }
138  else if (mxGetClassID(prhs[1]) == mxINT32_CLASS)
139  {
140  std::cout << "mxINT32_CLASS, int" << std::endl << std::flush;
141  callFunction<int>(plhs,prhs,nlhs,nrhs, isVectorImage);
142  }
143  else if (mxGetClassID(prhs[1]) == mxINT64_CLASS)
144  {
145  std::cout << "mxINT64_CLASS, int64_T" << std::endl << std::flush;
146  callFunction<int64_T>(plhs,prhs,nlhs,nrhs, isVectorImage);
147  }
148  else if (mxGetClassID(prhs[1]) == mxDOUBLE_CLASS)
149  {
150  std::cout << "mxDOUBLE_CLASS, double" << std::endl << std::flush;
151  callFunction<double>(plhs,prhs,nlhs,nrhs, isVectorImage);
152  }
153  else
154  utlGlobalException(true, "wrong type!!");
155 }
156 
mxArray * GetArrayStruct(const mxArray *pr_struct, const char *name)
Definition: mexutils.h:333
void SetITKImageInformation(itk::SmartPointer< ImageType > &image, const mxArray *originArray, const mxArray *spacingArray)
bool SaveImage(const SmartPointer< ImageType > &image, const std::string &filename, const std::string &printInfo="Writing Image:")
Definition: utlITK.h:157
void GetITKVectorImageFromMXArray(const mxArray *pr, SmartPointer< VectorImage< T, 3 > > &image)
Definition: mexITK.h:128
int mwSize
Definition: utils.h:38
const T & min(const T &a, const T &b)
Return the minimum between a and b.
Definition: utlCore.h:257
#define utlGlobalException(cond, expout)
Definition: utlCoreMacro.h:372
void callFunction(mxArray *plhs[], const mxArray *prhs[], const int nlhs, const int nrhs, const bool isVectorImage)
void GetString(const mxArray *pr, std::string &str)
Definition: mexutils.h:275
void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
void GetITKImageFromMXArray(const mxArray *pr, SmartPointer< Image< T, VImageDimension > > &image)
Definition: mexITK.h:55