DGtal 1.3.0
Loading...
Searching...
No Matches
ITKDicomReader.ih
1/**
2 * This program is free software: you can redistribute it and/or modify
3 * it under the terms of the GNU Lesser General Public License as
4 * published by the Free Software Foundation, either version 3 of the
5 * License, or (at your option) any later version.
6 *
7 * This program is distributed in the hope that it will be useful,
8 * but WITHOUT ANY WARRANTY; without even the implied warranty of
9 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
10 * GNU General Public License for more details.
11 *
12 * You should have received a copy of the GNU General Public License
13 * along with this program. If not, see <http://www.gnu.org/licenses/>.
14 *
15 **/
16
17/**
18 * @file ITKDicomReader.ih
19 * @author Boris Mansencal (\c boris.mansencal@labri.fr )
20 * LaBRI (CNRS, UMR 5800, University of Bordeaux, Bordeaux-INP, France
21 *
22 * @date 2019/02/05
23 *
24 * Header file for module ITKDicomReader.cpp
25 *
26 * This file is part of the DGtal library.
27 */
28
29#include "DGtal/images/ConstImageAdapter.h"
30#include "DGtal/io/readers/ITKReader.h"
31#if defined(__GNUG__)
32#pragma GCC diagnostic push
33#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
34#endif
35#if defined(__clang__)
36#pragma clang diagnostic push
37#pragma clang diagnostic ignored "-Wdocumentation"
38#endif
39#include <itkImageSeriesReader.h>
40#include <itkGDCMImageIO.h>
41#if defined(__clang__)
42#pragma clang diagnostic pop
43#endif
44#if defined(__GNUG__)
45#pragma GCC diagnostic pop
46#endif
47
48namespace DGtal {
49
50 template <typename I>
51 template <typename TFunctor>
52 I ITKDicomReader<I>::importDICOM( const std::vector<std::string> & filenames,
53 const TFunctor & aFunctor )
54 {
55 if ( filenames.empty() )
56 {
57 trace.error() << "[ITKDicomReader] empty filenames vector passed.";
58 throw IOException();
59 }
60
61 typedef typename Image::Domain Domain;
62 typedef itk::ImageIOBase::IOComponentType IOComponentType;
63 BOOST_CONCEPT_ASSERT( (concepts::CUnaryFunctor<TFunctor, ValueOut, Value>));
64 const IOComponentType componentType =
65 ITKReader<I>::getITKComponentType( filenames[0] );
66 // We suppose all file have the same 'componentType'.
67
68 switch ( componentType )
69 {
70
71 default:
72 case itk::ImageIOBase::UNKNOWNCOMPONENTTYPE:
73 {
74 trace.error()
75 << "[ITKDicomReader] Unknown and unsupported component type!";
76 throw IOException();
77 }
78 case itk::ImageIOBase::UCHAR:
79 {
80 typedef ImageContainerByITKImage<Domain, unsigned char> DGtalITKImage;
81 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
82 }
83 case itk::ImageIOBase::CHAR:
84 {
85 typedef ImageContainerByITKImage<Domain, char> DGtalITKImage;
86 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
87 }
88 case itk::ImageIOBase::USHORT:
89 {
90 typedef ImageContainerByITKImage<Domain, unsigned short> DGtalITKImage;
91 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
92 }
93 case itk::ImageIOBase::SHORT:
94 {
95 typedef ImageContainerByITKImage<Domain, short> DGtalITKImage;
96 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
97 }
98 case itk::ImageIOBase::UINT:
99 {
100 typedef ImageContainerByITKImage<Domain, unsigned int> DGtalITKImage;
101 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
102 }
103 case itk::ImageIOBase::INT:
104 {
105 typedef ImageContainerByITKImage<Domain, int> DGtalITKImage;
106 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
107 }
108 case itk::ImageIOBase::ULONG:
109 {
110 typedef ImageContainerByITKImage<Domain, unsigned long> DGtalITKImage;
111 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
112 }
113 case itk::ImageIOBase::LONG:
114 {
115 typedef ImageContainerByITKImage<Domain, long> DGtalITKImage;
116 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
117 }
118#if (ITK_VERSION_MAJOR > 4)\
119 || (ITK_VERSION_MAJOR == 4 && ITK_VERSION_MINOR >= 13)
120 case itk::ImageIOBase::ULONGLONG:
121 {
122 typedef ImageContainerByITKImage<Domain, unsigned long long>
123 DGtalITKImage;
124 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
125 }
126 case itk::ImageIOBase::LONGLONG:
127 {
128 typedef ImageContainerByITKImage<Domain, long long> DGtalITKImage;
129 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
130 }
131#endif
132 case itk::ImageIOBase::FLOAT:
133 {
134 typedef ImageContainerByITKImage<Domain, float> DGtalITKImage;
135 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
136 }
137 case itk::ImageIOBase::DOUBLE:
138 {
139 typedef ImageContainerByITKImage<Domain, double> DGtalITKImage;
140 return readDGtalImageFromITKtypes<DGtalITKImage>( filenames, aFunctor );
141 }
142 }
143 }
144
145
146 template <typename I>
147 template <typename Domain, typename PixelType>
148 inline
149 ImageContainerByITKImage<Domain, PixelType>
150 ITKDicomReader<I>::
151 importDicomFiles(const std::vector<std::string> & filenames)
152 {
153 typedef ImageContainerByITKImage<Domain, PixelType> ImageContainer;
154
155 const unsigned int dimension = Domain::dimension;
156 typedef itk::Image<PixelType, dimension> ItkImage;
157 typedef itk::ImageSeriesReader<ItkImage> ItkReader;
158
159 //typedef itk::GDCMImageIO ItkImageIO;
160
161 typedef typename ImageContainer::ITKImagePointer ITKImagePointer;
162 ITKImagePointer itkImage = nullptr;
163
164 try
165 {
166 typename ItkReader::Pointer reader = ItkReader::New();
167 //ItkImageIO::Pointer dicomIO = ItkImageIO::New();
168 //reader->SetImageIO( dicomIO );
169
170 reader->SetFileNames( filenames );
171 reader->Update();
172
173 itkImage = reader->GetOutput();
174 }
175 catch ( itk::ExceptionObject & e )
176 {
177 trace.error() << e;
178 throw IOException();
179 }
180 catch( ... )
181 {
182 trace.error() << "ITKDicomReader: can't read " << filenames.size()
183 << " files"<<std::endl;
184 throw IOException();
185 }
186
187 const typename ItkImage::SizeType& inputSize =
188 itkImage->GetLargestPossibleRegion().GetSize();
189 const unsigned int width = inputSize[0];
190 const unsigned int height = inputSize[1];
191 const unsigned int depth = inputSize[2];
192 if ( !height || !width || !depth )
193 {
194 trace.error() << "ITKDicomReader: one dimension is null (w=" << width
195 << ", h=" << height << ", d=" << depth << ")" << std::endl;
196 throw IOException();
197 }
198
199 const ImageContainer image( itkImage );
200
201 return image;
202 }
203
204
205 template <typename I>
206 template <typename Image, typename Domain, typename OrigValue,
207 typename TFunctor, typename Value>
208 Image
209 ITKDicomReader<I>::Aux<Image, Domain, OrigValue, TFunctor, Value>::
210 readDGtalImageFromITKtypes( const std::vector<std::string> & filenames,
211 const TFunctor & aFunctor )
212 {
213 typedef ImageContainerByITKImage<Domain, OrigValue> TypeDGtalImage;
214 TypeDGtalImage dgtalItkImage =
215 importDicomFiles<Domain, OrigValue>( filenames );
216
217 const Domain& domain = dgtalItkImage.domain();
218
219 typedef ConstImageAdapter<TypeDGtalImage, Domain, functors::Identity,
220 Value, TFunctor>
221 AdaptedImage;
222 const functors::Identity identityFunctor{};
223 const AdaptedImage adapted( dgtalItkImage, domain, identityFunctor,
224 aFunctor);
225
226 Image image( domain );
227 std::copy( adapted.constRange().begin(), adapted.constRange().end(),
228 image.range().outputIterator() );
229
230 return image;
231 }
232
233
234
235 //specialization
236 template <typename I>
237 template <typename Domain, typename OrigValue, typename TFunctor,
238 typename Value>
239 ImageContainerByITKImage<Domain, Value>
240 ITKDicomReader<I>::Aux<ImageContainerByITKImage<Domain, Value>, Domain,
241 OrigValue, TFunctor, Value>::
242 readDGtalImageFromITKtypes( const std::vector<std::string> & filenames,
243 const TFunctor & aFunctor )
244 {
245 typedef ImageContainerByITKImage<Domain, Value> Image;
246
247 typedef ImageContainerByITKImage<Domain, OrigValue> TypeDGtalImage;
248 TypeDGtalImage dgtalItkImage =
249 importDicomFiles<Domain, OrigValue>( filenames );
250
251 const Domain& domain = dgtalItkImage.domain();
252
253 typedef ConstImageAdapter<TypeDGtalImage, Domain, functors::Identity,
254 Value, TFunctor>
255 AdaptedImage;
256 const functors::Identity identityFunctor{};
257 const AdaptedImage adapted( dgtalItkImage, domain, identityFunctor,
258 aFunctor);
259
260 Image image( domain );
261 std::copy( adapted.constRange().begin(), adapted.constRange().end(),
262 image.range().outputIterator() );
263
264 //copy ITKImage spatial parameters
265 image.getITKImagePointer()->SetOrigin(
266 dgtalItkImage.getITKImagePointer()->GetOrigin() );
267 image.getITKImagePointer()->SetSpacing(
268 dgtalItkImage.getITKImagePointer()->GetSpacing() );
269 image.getITKImagePointer()->SetDirection(
270 dgtalItkImage.getITKImagePointer()->GetDirection() );
271
272 return image;
273 }
274
275
276 template <typename I>
277 template <typename TypeDGtalImage, typename TFunctor>
278 typename ITKDicomReader<I>::Image ITKDicomReader<I>::
279 readDGtalImageFromITKtypes( const std::vector<std::string> & filenames,
280 const TFunctor & aFunctor )
281 {
282 typedef typename Image::Domain Domain;
283 typedef typename TypeDGtalImage::Value OrigValue;
284
285 return Aux<Image, Domain, OrigValue, TFunctor, Value>::
286 readDGtalImageFromITKtypes( filenames, aFunctor );
287 }
288
289}//namespace