DGtal 1.3.0
Loading...
Searching...
No Matches
ITKReader.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
19 * @author Pierre Gueth (\c pierre.gueth@gmail.com )
20 * Laboratoire d'InfoRmatique en Image et Systèmes d'information - LIRIS (CNRS,
21 * UMR 5205), CNRS, France
22 *
23 *
24 * @author Bertrand Kerautret (\c bertrand.kerautret@loria.fr )
25 * LORIA (CNRS, UMR 7503), University of Lorraine, France
26 *
27 * @date 2013/10/28
28 *
29 * Header file for module ITKReader.cpp
30 *
31 * This file is part of the DGtal library.
32 */
33
34#include "DGtal/images/ConstImageAdapter.h"
35#if defined(__GNUG__)
36#pragma GCC diagnostic push
37#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
38#endif
39#if defined(__clang__)
40#pragma clang diagnostic push
41#pragma clang diagnostic ignored "-Wdocumentation"
42#endif
43#include <itkImageFileReader.h>
44#if defined(__clang__)
45#pragma clang diagnostic pop
46#endif
47#if defined(__GNUG__)
48#pragma GCC diagnostic pop
49#endif
50
51namespace DGtal {
52
53 template <typename I>
54 I ITKReader<I>::importITK(
55 const std::string & filename, bool shiftDomainUsingOrigin )
56 {
57 typename DGtal::ITKIOTrait<typename I::Value>::DefaultReadFunctor def;
58 return DGtal::ITKReader<I>::importITK(filename,def , shiftDomainUsingOrigin);
59 }
60
61 template <typename I>
62 template <typename TFunctor>
63 I ITKReader<I>::importITK(
64 const std::string & filename,
65 const TFunctor & aFunctor, bool shiftDomainUsingOrigin )
66 {
67 typedef typename Image::Domain Domain;
68 typedef itk::ImageIOBase::IOComponentType IOComponentType;
69 BOOST_CONCEPT_ASSERT( (concepts::CUnaryFunctor<TFunctor, ValueOut, Value>));
70 const IOComponentType componentType = getITKComponentType( filename );
71 switch ( componentType )
72 {
73
74 default:
75 case itk::ImageIOBase::UNKNOWNCOMPONENTTYPE:
76 trace.warning() << "[ITKReader] Unknown and unsupported component type! "
77 "File will be loaded with UCHAR component type" << std::endl;
78 //fallthrough
79 case itk::ImageIOBase::UCHAR:
80 {
81 typedef ImageContainerByITKImage<Domain, unsigned char> DGtalITKImage;
82 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
83 }
84 case itk::ImageIOBase::CHAR:
85 {
86 typedef ImageContainerByITKImage<Domain, char> DGtalITKImage;
87 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
88 }
89 case itk::ImageIOBase::USHORT:
90 {
91 typedef ImageContainerByITKImage<Domain, unsigned short> DGtalITKImage;
92 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
93 }
94 case itk::ImageIOBase::SHORT:
95 {
96 typedef ImageContainerByITKImage<Domain, short> DGtalITKImage;
97 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
98 }
99 case itk::ImageIOBase::UINT:
100 {
101 typedef ImageContainerByITKImage<Domain, unsigned int> DGtalITKImage;
102 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
103 }
104 case itk::ImageIOBase::INT:
105 {
106 typedef ImageContainerByITKImage<Domain, int> DGtalITKImage;
107 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
108 }
109 case itk::ImageIOBase::ULONG:
110 {
111 typedef ImageContainerByITKImage<Domain, unsigned long> DGtalITKImage;
112 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
113 }
114 case itk::ImageIOBase::LONG:
115 {
116 typedef ImageContainerByITKImage<Domain, long> DGtalITKImage;
117 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
118 }
119#if (ITK_VERSION_MAJOR > 4)\
120 || (ITK_VERSION_MAJOR == 4 && ITK_VERSION_MINOR >= 13)
121 case itk::ImageIOBase::ULONGLONG:
122 {
123 typedef ImageContainerByITKImage<Domain, unsigned long long>
124 DGtalITKImage;
125 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
126 }
127 case itk::ImageIOBase::LONGLONG:
128 {
129 typedef ImageContainerByITKImage<Domain, long long> DGtalITKImage;
130 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
131 }
132#endif
133 case itk::ImageIOBase::FLOAT:
134 {
135 typedef ImageContainerByITKImage<Domain, float> DGtalITKImage;
136 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
137 }
138 case itk::ImageIOBase::DOUBLE:
139 {
140 typedef ImageContainerByITKImage<Domain, double> DGtalITKImage;
141 return readDGtalImageFromITKtypes<DGtalITKImage>( filename, aFunctor, shiftDomainUsingOrigin );
142 }
143 }
144 }
145
146 template <typename I>
147 itk::ImageIOBase::IOComponentType
148 ITKReader<I>::getITKComponentType( const std::string & filename )
149 {
150 typedef itk::ImageIOBase::IOComponentType IOComponentType;
151 IOComponentType componentType;
152 try
153 {
154 itk::ImageIOBase::Pointer imageIO = itk::ImageIOFactory::CreateImageIO(
155 filename.c_str(), itk::ImageIOFactory::ReadMode );
156 if (imageIO)
157 {
158 imageIO->SetFileName( filename.c_str() );
159 imageIO->ReadImageInformation();
160 componentType = imageIO->GetComponentType();
161 }
162 }
163 catch ( itk::ExceptionObject & e )
164 {
165 trace.error() << e;
166 throw IOException();
167 }
168 return componentType;
169 }
170
171 template <typename I>
172 template <typename Domain, typename PixelType>
173 inline
174 ImageContainerByITKImage<Domain, PixelType>
175 ITKReader<I>::readDGtalITKImage(const std::string & filename, bool shiftDomainUsingOrigin)
176 {
177 typedef ImageContainerByITKImage<Domain, PixelType> TypeDGtalImage;
178 const unsigned int dimension = Domain::dimension;
179 typedef itk::Image<PixelType, dimension> ItkImage;
180 typedef itk::ImageFileReader<ItkImage> ITKImageReader;
181 typedef typename TypeDGtalImage::ITKImagePointer ITKImagePointer;
182 ITKImagePointer itk_image = nullptr;
183
184 try
185 {
186 typename ITKImageReader::Pointer reader = ITKImageReader::New();
187 reader->SetFileName( filename );
188 reader->Update();
189 reader->GetOutput();
190
191 itk_image = reader->GetOutput();
192 }
193 catch ( itk::ExceptionObject &e )
194 {
195 trace.error() << e;
196 throw IOException();
197 }
198
199 TypeDGtalImage image( itk_image );
200 if (shiftDomainUsingOrigin)
201 {
202 const auto c = itk_image->GetOrigin();
203 typename Domain::Point p;
204 for (unsigned int k = 0; k < dimension; k++)
205 {
206 p[k] = static_cast<typename Domain::Integer>(c[k]);
207 }
208 image.updateDomain(p);
209 }
210 return image;
211 }
212
213 template <typename I>
214 template <typename Image, typename Domain, typename OrigValue,
215 typename TFunctor, typename Value>
216 Image
217 ITKReader<I>::Aux<Image, Domain, OrigValue, TFunctor, Value>::
218 readDGtalImageFromITKtypes(const std::string & filename,
219 const TFunctor & aFunctor, bool shiftDomainUsingOrigin)
220 {
221 typedef ImageContainerByITKImage<Domain, OrigValue> TypeDGtalImage;
222 TypeDGtalImage dgtal_itk_image =
223 ITKReader<I>::readDGtalITKImage<Domain, OrigValue>( filename , shiftDomainUsingOrigin);
224
225 const Domain& domain = dgtal_itk_image.domain();
226
227 typedef ConstImageAdapter<TypeDGtalImage, Domain, functors::Identity,
228 Value, TFunctor> AdaptedImage;
229 const functors::Identity identityFunctor{};
230 const AdaptedImage adapted( dgtal_itk_image, domain, identityFunctor,
231 aFunctor );
232
233 Image image( domain );
234 std::copy( adapted.constRange().begin(), adapted.constRange().end(),
235 image.range().outputIterator() );
236
237 return image;
238 }
239
240
241 //specialization
242 template <typename I>
243 template <typename Domain, typename OrigValue, typename TFunctor,
244 typename Value>
245 ImageContainerByITKImage<Domain, Value>
246 ITKReader<I>::Aux<ImageContainerByITKImage<Domain, Value>, Domain, OrigValue,
247 TFunctor, Value>::
248 readDGtalImageFromITKtypes( const std::string & filename,
249 const TFunctor & aFunctor, bool shiftDomainUsingOrigin )
250 {
251 typedef ImageContainerByITKImage<Domain, Value> Image;
252
253 typedef ImageContainerByITKImage<Domain, OrigValue> TypeDGtalImage;
254 TypeDGtalImage dgtal_itk_image =
255 ITKReader<I>::readDGtalITKImage<Domain, OrigValue>( filename, shiftDomainUsingOrigin );
256
257 const Domain& domain = dgtal_itk_image.domain();
258
259 typedef ConstImageAdapter<TypeDGtalImage, Domain, functors::Identity,
260 Value, TFunctor>
261 AdaptedImage;
262 const functors::Identity identityFunctor{};
263 const AdaptedImage adapted( dgtal_itk_image, domain, identityFunctor,
264 aFunctor );
265
266 Image image( domain );
267 std::copy( adapted.constRange().begin(), adapted.constRange().end(),
268 image.range().outputIterator() );
269
270 //copy ITKImage spatial parameters
271 image.getITKImagePointer()->SetOrigin(
272 dgtal_itk_image.getITKImagePointer()->GetOrigin() );
273 image.getITKImagePointer()->SetSpacing(
274 dgtal_itk_image.getITKImagePointer()->GetSpacing() );
275 image.getITKImagePointer()->SetDirection(
276 dgtal_itk_image.getITKImagePointer()->GetDirection() );
277
278 return image;
279 }
280
281
282
283 template <typename I>
284 template <typename TypeDGtalImage, typename TFunctor>
285 typename ITKReader<I>::Image ITKReader<I>::
286 readDGtalImageFromITKtypes( const std::string & filename,
287 const TFunctor & aFunctor, bool shiftDomainUsingOrigin )
288 {
289 typedef typename Image::Domain Domain;
290 typedef typename TypeDGtalImage::Value OrigValue;
291
292 return Aux<Image, Domain, OrigValue, TFunctor,
293 Value>::readDGtalImageFromITKtypes( filename, aFunctor, shiftDomainUsingOrigin );
294 }
295
296}//namespace