ITK  6.0.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  * https://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 
64 
66  template <typename TOutputImageType>
67  static typename TOutputImageType::Pointer
68  IplImageToITKImage(const IplImage * in);
69 
71  template <typename TOutputImageType>
72  static typename TOutputImageType::Pointer
73  CVMatToITKImage(const cv::Mat & in);
74 
76  template <typename TInputImageType>
77  static IplImage *
78  ITKImageToIplImage(const TInputImageType * in, bool force3Channels = false);
79 
81  template <typename TInputImageType>
82  static cv::Mat
83  ITKImageToCVMat(const TInputImageType * in, bool force3Channels = false);
84 
85 private:
92  template <typename TOutputImageType, typename TPixel>
93  static void
94  ITKConvertIplImageBuffer(const IplImage * in, TOutputImageType * out, unsigned int iDepth)
95  {
96  // Typedefs
97  using ImageType = TOutputImageType;
98  using OutputPixelType = typename ImageType::PixelType;
101  unsigned int inChannels = in->nChannels;
103 
104  // Manage unsupported conversions.
105  if ((inChannels != outChannels || (inChannels == 3 && outChannels == 3)) &&
106  (iDepth == IPL_DEPTH_8S || iDepth == IPL_DEPTH_16S || iDepth == IPL_DEPTH_32S || iDepth == IPL_DEPTH_64F))
107  {
108  itkGenericExceptionMacro("OpenCV IplImage to ITK Image conversion - the necessary color"
109  " conversion is not supported for the input OpenCV pixel type");
110  }
111 
112  // Manage input/output types mismatch
113  checkMatchingTypes<TPixel, OutputPixelType>(outChannels);
114 
115  // We only change current if it no longer points at in, so this is safe
116  auto * current = const_cast<IplImage *>(in);
117 
118  bool freeCurrent = false;
119  if (inChannels == 3 && outChannels == 1)
120  {
121  current = cvCreateImage(cvSize(in->width, in->height), iDepth, 1);
122  cvCvtColor(in, current, CV_BGR2GRAY);
123  freeCurrent = true;
124  }
125  else if (inChannels == 1 && outChannels == 3)
126  {
127  current = cvCreateImage(cvSize(in->width, in->height), iDepth, 3);
128  cvCvtColor(in, current, CV_GRAY2RGB);
129  freeCurrent = true;
130  }
131  else if (inChannels == 3 && outChannels == 3)
132  {
133  current = cvCreateImage(cvSize(in->width, in->height), iDepth, 3);
134  cvCvtColor(in, current, CV_BGR2RGB);
135  freeCurrent = true;
136  }
137  else if (inChannels != 1 || outChannels != 1)
138  {
139  itkGenericExceptionMacro("Conversion from " << inChannels << " channels to " << outChannels
140  << " channels is not supported");
141  }
142 
143  ITKConvertImageBuffer<TOutputImageType, TPixel>(
144  current->imageData, out, current->nChannels, current->width, current->height, current->widthStep);
145 
146  if (freeCurrent)
147  {
148  cvReleaseImage(&current);
149  }
150  }
151 
152  template <typename TOutputImageType, typename TPixel>
153  static void
154  ITKConvertMatImageBuffer(const cv::Mat & in, TOutputImageType * out)
155  {
156  using namespace cv;
157 
158  // Typedefs
159  using ImageType = TOutputImageType;
160  using OutputPixelType = typename ImageType::PixelType;
161 
162  unsigned int inChannels = in.channels();
163  unsigned int iDepth = in.depth();
165 
166  // Manage pixel/component types mismatch and impossible conversions.
167  if ((inChannels != outChannels || (inChannels == 3 && outChannels == 3)) &&
168  (iDepth == CV_8S || iDepth == CV_16S || iDepth == CV_32S || iDepth == CV_64F))
169  {
170  itkGenericExceptionMacro("OpenCV Mat to ITK Image conversion - the necessary color"
171  " conversion is not supported for the input OpenCV pixel type");
172  }
173 
174  // Manage input/output types mismatch
175  checkMatchingTypes<TPixel, OutputPixelType>(outChannels);
176 
177  Mat current;
178 
179  if (inChannels == 3 && outChannels == 1)
180  {
181  cvtColor(in, current, COLOR_BGR2GRAY);
182  }
183  else if (inChannels == 1 && outChannels == 3)
184  {
185  cvtColor(in, current, COLOR_GRAY2RGB);
186  }
187  else if (inChannels == 3 && outChannels == 3)
188  {
189  cvtColor(in, current, COLOR_BGR2RGB);
190  }
191  else if (inChannels != 1 || outChannels != 1)
192  {
193  itkGenericExceptionMacro("Conversion from " << inChannels << " channels to " << outChannels
194  << " channels is not supported");
195  }
196  else
197  {
198  current = in;
199  }
200 
201  ITKConvertImageBuffer<TOutputImageType, TPixel>(
202  reinterpret_cast<char *>(current.ptr()), out, current.channels(), current.cols, current.rows, current.step);
203  }
204 
205  template <typename InputPixelType, typename OutputPixelType>
206  static void
207  checkMatchingTypes(unsigned int outChannels)
208  {
209  if (outChannels == 3)
210  {
211  if (typeid(RGBPixel<InputPixelType>) != typeid(OutputPixelType))
212  {
213  itkGenericExceptionMacro("OpenCV to ITK conversion channel component type mismatch");
214  }
215  }
216  else if (typeid(InputPixelType) != typeid(OutputPixelType))
217  {
218  itkGenericExceptionMacro("OpenCV to ITK conversion pixel type mismatch");
219  }
220  }
221 
222  template <typename TOutputImageType, typename TPixel>
223  static void
224  ITKConvertImageBuffer(const char * in,
225  TOutputImageType * out,
226  unsigned int inChannels,
227  int imgWidth,
228  int imgHeight,
229  int widthStep)
230  {
231  using ImageType = TOutputImageType;
232  using OutputPixelType = typename ImageType::PixelType;
233  using ConvertPixelTraits = DefaultConvertPixelTraits<OutputPixelType>;
234 
235  bool isVectorImage(strcmp(out->GetNameOfClass(), "VectorImage") == 0);
236 
237  typename ImageType::RegionType region;
238  typename ImageType::RegionType::SizeType size;
239  typename ImageType::RegionType::IndexType start;
240  typename ImageType::SpacingType spacing;
241  size.Fill(1);
242  size[0] = imgWidth;
243  size[1] = imgHeight;
244  start.Fill(0);
245  spacing.Fill(1);
246  region.SetSize(size);
247  region.SetIndex(start);
248  out->SetRegions(region);
249  out->SetSpacing(spacing);
250  out->Allocate();
251  size_t lineLength = imgWidth * inChannels * sizeof(TPixel);
252  auto * unpaddedBuffer = reinterpret_cast<void *>(new TPixel[imgHeight * lineLength]);
253  unsigned int paddedBufPos = 0;
254  unsigned int unpaddedBufPos = 0;
255 
256  for (int i = 0; i < imgHeight; ++i)
257  {
258  memcpy(&(reinterpret_cast<char *>(unpaddedBuffer)[unpaddedBufPos]), in + paddedBufPos, lineLength);
259  paddedBufPos += widthStep;
260  unpaddedBufPos += lineLength;
261  }
262 
263  if (isVectorImage)
264  {
266  static_cast<TPixel *>(unpaddedBuffer),
267  inChannels,
268  out->GetPixelContainer()->GetBufferPointer(),
269  out->GetPixelContainer()->Size());
270  }
271  else
272  {
274  static_cast<TPixel *>(unpaddedBuffer),
275  inChannels,
276  out->GetPixelContainer()->GetBufferPointer(),
277  out->GetPixelContainer()->Size());
278  }
279 
280  delete[] reinterpret_cast<TPixel *>(unpaddedBuffer);
281  }
282 
283  template <typename TPixel, unsigned int VDimension>
285  {
286  static void
287  Padding(const Image<TPixel, VDimension> * itkNotUsed(in), IplImage * itkNotUsed(out))
288  {}
289  };
290 
291  template <typename TValue, unsigned int VDimension>
292  struct HandleRGBPixel<RGBPixel<TValue>, VDimension>
293  {
294  using ValueType = TValue;
297 
298  static void
299  Padding(const ImageType * in, IplImage * out)
300  {
301  typename ImageType::IndexType pixelIndex = { { 0, 0 } };
302 
303  for (int r = 0; r < out->height; ++r)
304  {
305  auto * ptr = reinterpret_cast<ValueType *>(out->imageData + r * out->widthStep);
306  for (int c = 0; c < out->width; ++c)
307  {
308  pixelIndex[0] = c;
309  pixelIndex[1] = r;
310  typename ImageType::PixelType pixel = in->GetPixel(pixelIndex);
311 
312  for (unsigned int i = 0; i < 3; ++i)
313  {
314  *ptr++ = pixel[i];
315  }
316  }
317  }
318  }
319  };
320 }; // end class OpenCVImageBridge
321 
322 } // end namespace itk
323 
324 #ifndef ITK_MANUAL_INSTANTIATION
325 # include "itkOpenCVImageBridge.hxx"
326 #endif
327 
328 #endif
itk::OpenCVImageBridge
This class provides static methods to convert between OpenCV images and itk::Image.
Definition: itkOpenCVImageBridge.h:58
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Index< VImageDimension >
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:294
itkConvertPixelBuffer.h
itk::OpenCVImageBridge::HandleRGBPixel::Padding
static void Padding(const Image< TPixel, VDimension > *, IplImage *)
Definition: itkOpenCVImageBridge.h:287
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:219
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::OpenCVImageBridge::checkMatchingTypes
static void checkMatchingTypes(unsigned int outChannels)
Definition: itkOpenCVImageBridge.h:207
itkImage.h
itk::OpenCVImageBridge::ITKImageToIplImage
static IplImage * ITKImageToIplImage(const TInputImageType *in, bool force3Channels=false)
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:272
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::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:154
itk::OpenCVImageBridge::ITKConvertImageBuffer
static void ITKConvertImageBuffer(const char *in, TOutputImageType *out, unsigned int inChannels, int imgWidth, int imgHeight, int widthStep)
Definition: itkOpenCVImageBridge.h:224
itk::OpenCVImageBridge::ITKConvertIplImageBuffer
static void ITKConvertIplImageBuffer(const IplImage *in, TOutputImageType *out, unsigned int iDepth)
Definition: itkOpenCVImageBridge.h:94
itk::Image::PixelType
TPixel PixelType
Definition: itkImage.h:108
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itk::ConvertPixelBuffer::ConvertVectorImage
static void ConvertVectorImage(const 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:299
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::OpenCVImageBridge::CVMatToITKImage
static TOutputImageType::Pointer CVMatToITKImage(const cv::Mat &in)
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::ConvertPixelBuffer::Convert
static void Convert(const InputPixelType *inputData, int inputNumberOfComponents, OutputPixelType *outputData, vcl_size_t size)
itk::OpenCVImageBridge::HandleRGBPixel
Definition: itkOpenCVImageBridge.h:284