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