ITK  5.2.0
Insight Toolkit
itkOpenCVImageBridge.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkOpenCVImageBridge_h
19 #define itkOpenCVImageBridge_h
20 
21 #include <string>
22 
23 #include "itkImage.h"
25 #include "itkConvertPixelBuffer.h"
26 
27 #include "opencv2/core/version.hpp"
28 #if !defined(CV_VERSION_EPOCH)
29 // OpenCV 3+
30 # include "opencv2/imgproc.hpp"
31 # include "opencv2/imgproc/types_c.h" // CV_RGB2BGR, CV_BGR2GRAY, ...
32 # include "opencv2/imgproc/imgproc_c.h" // cvCvtColor
33 #else
34 // OpenCV 2.4.x
35 # include "cv.h"
36 # include "highgui.h"
37 #endif
38 
39 namespace itk
40 {
41 
59 {
60 public:
61  ITK_DISALLOW_COPY_AND_MOVE(OpenCVImageBridge);
62 
65 
67  template <typename TOutputImageType>
68  static typename TOutputImageType::Pointer
69  IplImageToITKImage(const IplImage * in);
70 
72  template <typename TOutputImageType>
73  static typename TOutputImageType::Pointer
74  CVMatToITKImage(const cv::Mat & in);
75 
77  template <typename TInputImageType>
78  static IplImage *
79  ITKImageToIplImage(const TInputImageType * in, bool force3Channels = false);
80 
82  template <typename TInputImageType>
83  static cv::Mat
84  ITKImageToCVMat(const TInputImageType * in, bool force3Channels = false);
85 
86 private:
93  template <typename TOutputImageType, typename TPixel>
94  static void
95  ITKConvertIplImageBuffer(const IplImage * in, TOutputImageType * out, unsigned int iDepth)
96  {
97  // Typedefs
98  using ImageType = TOutputImageType;
99  using OutputPixelType = typename ImageType::PixelType;
101 
102  unsigned int inChannels = in->nChannels;
104 
105  // Manage unsupported conversions.
106  if ((inChannels != outChannels || (inChannels == 3 && outChannels == 3)) &&
107  (iDepth == IPL_DEPTH_8S || iDepth == IPL_DEPTH_16S || iDepth == IPL_DEPTH_32S || iDepth == IPL_DEPTH_64F))
108  {
109  itkGenericExceptionMacro("OpenCV IplImage to ITK Image conversion - the necessary color"
110  " conversion is not supported for the input OpenCV pixel type");
111  }
112 
113  // Manage input/output types mismatch
114  checkMatchingTypes<TPixel, OutputPixelType>(outChannels);
115 
116  // We only change current if it no longer points at in, so this is safe
117  IplImage * current = const_cast<IplImage *>(in);
118 
119  bool freeCurrent = false;
120  if (inChannels == 3 && outChannels == 1)
121  {
122  current = cvCreateImage(cvSize(in->width, in->height), iDepth, 1);
123  cvCvtColor(in, current, CV_BGR2GRAY);
124  freeCurrent = true;
125  }
126  else if (inChannels == 1 && outChannels == 3)
127  {
128  current = cvCreateImage(cvSize(in->width, in->height), iDepth, 3);
129  cvCvtColor(in, current, CV_GRAY2RGB);
130  freeCurrent = true;
131  }
132  else if (inChannels == 3 && outChannels == 3)
133  {
134  current = cvCreateImage(cvSize(in->width, in->height), iDepth, 3);
135  cvCvtColor(in, current, CV_BGR2RGB);
136  freeCurrent = true;
137  }
138  else if (inChannels != 1 || outChannels != 1)
139  {
140  itkGenericExceptionMacro("Conversion from " << inChannels << " channels to " << outChannels
141  << " channels is not supported");
142  }
143 
144  ITKConvertImageBuffer<TOutputImageType, TPixel>(
145  current->imageData, out, current->nChannels, current->width, current->height, current->widthStep);
146 
147  if (freeCurrent)
148  {
149  cvReleaseImage(&current);
150  }
151  }
152 
153  template <typename TOutputImageType, typename TPixel>
154  static void
155  ITKConvertMatImageBuffer(const cv::Mat & in, TOutputImageType * out)
156  {
157  using namespace cv;
158 
159  // Typedefs
160  using ImageType = TOutputImageType;
161  using OutputPixelType = typename ImageType::PixelType;
162 
163  unsigned int inChannels = in.channels();
164  unsigned int iDepth = in.depth();
166 
167  // Manage pixel/component types mismatch and impossible conversions.
168  if ((inChannels != outChannels || (inChannels == 3 && outChannels == 3)) &&
169  (iDepth == CV_8S || iDepth == CV_16S || iDepth == CV_32S || iDepth == CV_64F))
170  {
171  itkGenericExceptionMacro("OpenCV Mat to ITK Image conversion - the necessary color"
172  " conversion is not supported for the input OpenCV pixel type");
173  }
174 
175  // Manage input/output types mismatch
176  checkMatchingTypes<TPixel, OutputPixelType>(outChannels);
177 
178  Mat current;
179 
180  if (inChannels == 3 && outChannels == 1)
181  {
182  cvtColor(in, current, COLOR_BGR2GRAY);
183  }
184  else if (inChannels == 1 && outChannels == 3)
185  {
186  cvtColor(in, current, COLOR_GRAY2RGB);
187  }
188  else if (inChannels == 3 && outChannels == 3)
189  {
190  cvtColor(in, current, COLOR_BGR2RGB);
191  }
192  else if (inChannels != 1 || outChannels != 1)
193  {
194  itkGenericExceptionMacro("Conversion from " << inChannels << " channels to " << outChannels
195  << " channels is not supported");
196  }
197  else
198  {
199  current = in;
200  }
201 
202  ITKConvertImageBuffer<TOutputImageType, TPixel>(
203  reinterpret_cast<char *>(current.ptr()), out, current.channels(), current.cols, current.rows, current.step);
204  }
205 
206  template <typename InputPixelType, typename OutputPixelType>
207  static void
208  checkMatchingTypes(unsigned int outChannels)
209  {
210  if (outChannels == 3)
211  {
212  if (typeid(RGBPixel<InputPixelType>) != typeid(OutputPixelType))
213  {
214  itkGenericExceptionMacro("OpenCV to ITK conversion channel component type mismatch");
215  }
216  }
217  else if (typeid(InputPixelType) != typeid(OutputPixelType))
218  {
219  itkGenericExceptionMacro("OpenCV to ITK conversion pixel type mismatch");
220  }
221  }
222 
223  template <typename TOutputImageType, typename TPixel>
224  static void
225  ITKConvertImageBuffer(const char * in,
226  TOutputImageType * out,
227  unsigned int inChannels,
228  int imgWidth,
229  int imgHeight,
230  int widthStep)
231  {
232  using ImageType = TOutputImageType;
233  using OutputPixelType = typename ImageType::PixelType;
234  using ConvertPixelTraits = DefaultConvertPixelTraits<OutputPixelType>;
235 
236  bool isVectorImage(strcmp(out->GetNameOfClass(), "VectorImage") == 0);
237 
238  typename ImageType::RegionType region;
239  typename ImageType::RegionType::SizeType size;
240  typename ImageType::RegionType::IndexType start;
241  typename ImageType::SpacingType spacing;
242  size.Fill(1);
243  size[0] = imgWidth;
244  size[1] = imgHeight;
245  start.Fill(0);
246  spacing.Fill(1);
247  region.SetSize(size);
248  region.SetIndex(start);
249  out->SetRegions(region);
250  out->SetSpacing(spacing);
251  out->Allocate();
252  size_t lineLength = imgWidth * inChannels * sizeof(TPixel);
253  void * unpaddedBuffer = reinterpret_cast<void *>(new TPixel[imgHeight * lineLength]);
254  unsigned int paddedBufPos = 0;
255  unsigned int unpaddedBufPos = 0;
256 
257  for (int i = 0; i < imgHeight; ++i)
258  {
259  memcpy(&(reinterpret_cast<char *>(unpaddedBuffer)[unpaddedBufPos]), in + paddedBufPos, lineLength);
260  paddedBufPos += widthStep;
261  unpaddedBufPos += lineLength;
262  }
263 
264  if (isVectorImage)
265  {
267  static_cast<TPixel *>(unpaddedBuffer),
268  inChannels,
269  out->GetPixelContainer()->GetBufferPointer(),
270  out->GetPixelContainer()->Size());
271  }
272  else
273  {
275  static_cast<TPixel *>(unpaddedBuffer),
276  inChannels,
277  out->GetPixelContainer()->GetBufferPointer(),
278  out->GetPixelContainer()->Size());
279  }
280 
281  delete[] reinterpret_cast<TPixel *>(unpaddedBuffer);
282  }
283 
284  template <typename TPixel, unsigned int VDimension>
286  {
287  static void
288  Padding(const Image<TPixel, VDimension> * itkNotUsed(in), IplImage * itkNotUsed(out))
289  {}
290  };
291 
292  template <typename TValue, unsigned int VDimension>
293  struct HandleRGBPixel<RGBPixel<TValue>, VDimension>
294  {
295  using ValueType = TValue;
298 
299  static void
300  Padding(const ImageType * in, IplImage * out)
301  {
302  typename ImageType::IndexType pixelIndex = { { 0, 0 } };
303 
304  for (int r = 0; r < out->height; r++)
305  {
306  ValueType * ptr = reinterpret_cast<ValueType *>(out->imageData + r * out->widthStep);
307  for (int c = 0; c < out->width; c++)
308  {
309  pixelIndex[0] = c;
310  pixelIndex[1] = r;
311  typename ImageType::PixelType pixel = in->GetPixel(pixelIndex);
312 
313  for (unsigned int i = 0; i < 3; i++)
314  {
315  *ptr++ = pixel[i];
316  }
317  }
318  }
319  }
320  };
321 }; // end class OpenCVImageBridge
322 
323 } // end namespace itk
324 
325 #ifndef ITK_MANUAL_INSTANTIATION
326 # include "itkOpenCVImageBridge.hxx"
327 #endif
328 
329 #endif
itk::OpenCVImageBridge
This class provides static methods to convert between OpenCV images and itk::Image.
Definition: itkOpenCVImageBridge.h:58
itk::RGBPixel
Represent Red, Green and Blue components for color images.
Definition: itkRGBPixel.h:58
itk::OpenCVImageBridge::HandleRGBPixel< RGBPixel< TValue >, VDimension >::ValueType
TValue ValueType
Definition: itkOpenCVImageBridge.h:295
itkConvertPixelBuffer.h
itk::OpenCVImageBridge::HandleRGBPixel::Padding
static void Padding(const Image< TPixel, VDimension > *, IplImage *)
Definition: itkOpenCVImageBridge.h:288
itk::OpenCVImageBridge::IplImageToITKImage
static TOutputImageType::Pointer IplImageToITKImage(const IplImage *in)
itk::Image::GetPixel
const TPixel & GetPixel(const IndexType &index) const
Get a pixel (read only version).
Definition: itkImage.h:217
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::OpenCVImageBridge::checkMatchingTypes
static void checkMatchingTypes(unsigned int outChannels)
Definition: itkOpenCVImageBridge.h:208
itkImage.h
itk::OpenCVImageBridge::ITKImageToIplImage
static IplImage * ITKImageToIplImage(const TInputImageType *in, bool force3Channels=false)
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:270
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itkDefaultConvertPixelTraits.h
itk::OpenCVImageBridge::ITKImageToCVMat
static cv::Mat ITKImageToCVMat(const TInputImageType *in, bool force3Channels=false)
itk::ConvertPixelBuffer::ConvertVectorImage
static void ConvertVectorImage(InputPixelType *inputData, int inputNumberOfComponents, OutputPixelType *outputData, vcl_size_t size)
itk::DefaultConvertPixelTraits
Traits class used to by ConvertPixels to convert blocks of pixels.
Definition: itkDefaultConvertPixelTraits.h:41
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::OpenCVImageBridge::ITKConvertMatImageBuffer
static void ITKConvertMatImageBuffer(const cv::Mat &in, TOutputImageType *out)
Definition: itkOpenCVImageBridge.h:155
itk::OpenCVImageBridge::ITKConvertImageBuffer
static void ITKConvertImageBuffer(const char *in, TOutputImageType *out, unsigned int inChannels, int imgWidth, int imgHeight, int widthStep)
Definition: itkOpenCVImageBridge.h:225
itk::OpenCVImageBridge::ITKConvertIplImageBuffer
static void ITKConvertIplImageBuffer(const IplImage *in, TOutputImageType *out, unsigned int iDepth)
Definition: itkOpenCVImageBridge.h:95
itk::Image::PixelType
TPixel PixelType
Definition: itkImage.h:106
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:58
itk::ConvertPixelBuffer::Convert
static void Convert(InputPixelType *inputData, int inputNumberOfComponents, OutputPixelType *outputData, vcl_size_t size)
itk::OpenCVImageBridge::HandleRGBPixel< RGBPixel< TValue >, VDimension >::Padding
static void Padding(const ImageType *in, IplImage *out)
Definition: itkOpenCVImageBridge.h:300
itk::Image::IndexType
typename Superclass::IndexType IndexType
Definition: itkImage.h:132
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::OpenCVImageBridge::CVMatToITKImage
static TOutputImageType::Pointer CVMatToITKImage(const cv::Mat &in)
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::OpenCVImageBridge::HandleRGBPixel
Definition: itkOpenCVImageBridge.h:285