ITK  5.4.0
Insight Toolkit
itkImageAdaptor.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 itkImageAdaptor_h
19#define itkImageAdaptor_h
20
21#include "itkImage.h"
22
23namespace itk
24{
25
26template <typename TPixelType, unsigned int VImageDimension>
27class VectorImage;
28
55template <typename TImage, typename TAccessor>
56class ITK_TEMPLATE_EXPORT ImageAdaptor : public ImageBase<TImage::ImageDimension>
57{
58public:
59 ITK_DISALLOW_COPY_AND_MOVE(ImageAdaptor);
60
65 static constexpr unsigned int ImageDimension = TImage::ImageDimension;
66
73
75 itkOverrideGetNameOfClassMacro(ImageAdaptor);
76
78 using InternalImageType = TImage;
79
81 itkNewMacro(Self);
82
85 using PixelType = typename TAccessor::ExternalType;
86
89 using InternalPixelType = typename TAccessor::InternalType;
90
92
95 using AccessorType = TAccessor;
96
99 using AccessorFunctorType = typename InternalImageType::AccessorFunctorType::template Rebind<Self>::Type;
100
102 using typename Superclass::IndexType;
104
106 using typename Superclass::SizeType;
108
110 using typename Superclass::OffsetType;
112
115 using typename Superclass::RegionType;
116
119 using typename Superclass::SpacingType;
120
123 using typename Superclass::PointType;
124
128 using typename Superclass::DirectionType;
129
130
137 template <typename UPixelType, unsigned int UImageDimension = TImage::ImageDimension>
138 struct Rebind
139 {
141 };
142
143 template <typename UPixelType, unsigned int VUImageDimension = TImage::ImageDimension>
145
146
153 void
154 SetLargestPossibleRegion(const RegionType & region) override;
155
159 void
160 SetBufferedRegion(const RegionType & region) override;
161
165 void
166 SetRequestedRegion(const RegionType & region) override;
167
172 void
173 SetRequestedRegion(const DataObject * data) override;
174
181 const RegionType &
182 GetRequestedRegion() const override;
183
192 const RegionType &
193 GetLargestPossibleRegion() const override;
194
200 const RegionType &
201 GetBufferedRegion() const override;
202
204 void
205 Allocate(bool initialize = false) override;
206
209 void
210 Initialize() override;
211
213 void
214 SetPixel(const IndexType & index, const PixelType & value)
215 {
216 m_PixelAccessor.Set(m_Image->GetPixel(index), value);
217 }
218
220 PixelType
221 GetPixel(const IndexType & index) const
222 {
223 return m_PixelAccessor.Get(m_Image->GetPixel(index));
224 }
225
227 PixelType operator[](const IndexType & index) const { return m_PixelAccessor.Get(m_Image->GetPixel(index)); }
228
230 const OffsetValueType *
232
236
239 using PixelContainer = typename TImage::PixelContainer;
240 using PixelContainerPointer = typename TImage::PixelContainerPointer;
241 using PixelContainerConstPointer = typename TImage::PixelContainerConstPointer;
242
246 {
247 return m_Image->GetPixelContainer();
248 }
249
250 const PixelContainer *
252 {
253 return m_Image->GetPixelContainer();
254 }
255
258 void
260
271 virtual void
272 Graft(const Self * imgData);
273
276
281
282 const InternalPixelType *
284
286 void
287 SetSpacing(const SpacingType & spacing) override;
288
289 void
290 SetSpacing(const double * spacing /*[ImageDimension]*/) override;
291
292 void
293 SetSpacing(const float * spacing /*[ImageDimension]*/) override;
294
298 const SpacingType &
299 GetSpacing() const override;
300
304 const PointType &
305 GetOrigin() const override;
306
308 void
309 SetOrigin(const PointType origin) override;
310
311 void
312 SetOrigin(const double * origin /*[ImageDimension]*/) override;
313
314 void
315 SetOrigin(const float * origin /*[ImageDimension]*/) override;
316
318 void
319 SetDirection(const DirectionType & direction) override;
320
324 const DirectionType &
325 GetDirection() const override;
326
328 virtual void
329 SetImage(TImage *);
330
332 void
333 Modified() const override;
334
337 GetMTime() const override;
338
342 {
343 return m_PixelAccessor;
344 }
345
347 const AccessorType &
349 {
350 return m_PixelAccessor;
351 }
352
354 void
356 {
357 m_PixelAccessor = accessor;
358 }
359
361 void
362 Update() override;
363
364 void
365 CopyInformation(const DataObject * data) override;
366
369 void
371
372 void
374
375 void
377
378 void
380
381 bool
383
385 template <typename TIndexRep, typename TCoordRep>
388 {
389 return m_Image->template TransformPhysicalPointToContinuousIndex<TIndexRep>(point);
390 }
391
400 template <typename TCoordRep>
401 ITK_NODISCARD("Call the overload which has the point as the only parameter and returns the index")
402 bool TransformPhysicalPointToContinuousIndex(const Point<TCoordRep, Self::ImageDimension> & point,
403 ContinuousIndex<TCoordRep, Self::ImageDimension> & index) const
404 {
405 return m_Image->TransformPhysicalPointToContinuousIndex(point, index);
406 }
407
409 template <typename TCoordRep>
410 [[nodiscard]] IndexType
412 {
413 return m_Image->TransformPhysicalPointToIndex(point);
414 }
415
424 template <typename TCoordRep>
425 ITK_NODISCARD("Call the overload which has the point as the only parameter and returns the index")
426 bool TransformPhysicalPointToIndex(const Point<TCoordRep, Self::ImageDimension> & point, IndexType & index) const
427 {
428 return m_Image->TransformPhysicalPointToIndex(point, index);
429 }
430
435 template <typename TCoordRep>
436 void
439 {
440 m_Image->TransformContinuousIndexToPhysicalPoint(index, point);
441 }
442
444 template <typename TCoordRep, typename TIndexRep>
447 {
448 return m_Image->template TransformContinuousIndexToPhysicalPoint<TIndexRep>(index);
449 }
450
456 template <typename TCoordRep>
457 void
459 {
460 m_Image->TransformIndexToPhysicalPoint(index, point);
461 }
462
464 template <typename TCoordRep>
467 {
468 return m_Image->template TransformIndexToPhysicalPoint<TCoordRep>(index);
469 }
470
471 template <typename TCoordRep>
472 void
474 FixedArray<TCoordRep, Self::ImageDimension> & outputGradient) const
475 {
476 m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
477 }
478
479 template <typename TVector>
480 [[nodiscard]] TVector
481 TransformLocalVectorToPhysicalVector(const TVector & inputGradient) const
482 {
483 TVector outputGradient;
484 TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
485 return outputGradient;
486 }
487
488 template <typename TCoordRep>
489 void
491 FixedArray<TCoordRep, Self::ImageDimension> & outputGradient) const
492 {
493 m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
494 }
495
496 template <typename TVector>
497 [[nodiscard]] TVector
498 TransformPhysicalVectorToLocalVector(const TVector & inputGradient) const
499 {
500 TVector outputGradient;
501 TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
502 return outputGradient;
503 }
504
505protected:
507 ~ImageAdaptor() override = default;
508 void
509 PrintSelf(std::ostream & os, Indent indent) const override;
510 void
511 Graft(const DataObject * data) override;
512 using Superclass::Graft;
513
514private:
515 // a specialized method to update PixelAccessors for VectorImages,
516 // to have the correct vector length of the image.
517 template <typename TPixelType>
518 void
520 {
521 this->m_PixelAccessor.SetVectorLength(this->m_Image->GetNumberOfComponentsPerPixel());
522 }
523
524 // The other image types don't expect an accessor which needs any updates
525 template <typename T>
526 void
527 UpdateAccessor(T * itkNotUsed(dummy))
528 {}
529
530 // Adapted image, most of the calls to ImageAdaptor
531 // will be delegated to this image
532 typename TImage::Pointer m_Image{};
533
534 // Data accessor object,
535 // it converts the presentation of a pixel
536 AccessorType m_PixelAccessor{};
537};
538} // end namespace itk
539
540#ifndef ITK_MANUAL_INSTANTIATION
541# include "itkImageAdaptor.hxx"
542#endif
543
544#endif
A templated class holding a point in n-Dimensional image space.
Base class for all data objects in ITK.
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:54
Give access to partial aspects of voxels from an Image.
void SetOrigin(const float *origin) override
void CopyInformation(const DataObject *data) override
typename InternalImageType::AccessorFunctorType::template Rebind< Self >::Type AccessorFunctorType
typename OffsetType::OffsetValueType OffsetValueType
void SetOrigin(const double *origin) override
void SetPixelContainer(PixelContainer *container)
const PointType & GetOrigin() const override
void Modified() const override
void SetPixel(const IndexType &index, const PixelType &value)
void Initialize() override
PixelContainerPointer GetPixelContainer()
void SetRequestedRegion(const DataObject *data) override
typename TImage::PixelContainerPointer PixelContainerPointer
InternalPixelType * InternalPixelPointerType
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, Self::ImageDimension > &inputGradient, FixedArray< TCoordRep, Self::ImageDimension > &outputGradient) const
typename SizeType::SizeValueType SizeValueType
void PropagateRequestedRegion() override
void UpdateAccessor(typename itk::VectorImage< TPixelType, ImageDimension > *)
void SetRequestedRegion(const RegionType &region) override
const PixelContainer * GetPixelContainer() const
PixelType operator[](const IndexType &index) const
Point< TCoordRep, TImage::ImageDimension > TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, Self::ImageDimension > &index) const
typename IndexType::IndexValueType IndexValueType
const DirectionType & GetDirection() const override
bool VerifyRequestedRegion() override
void UpdateOutputInformation() override
void SetSpacing(const SpacingType &spacing) override
void SetBufferedRegion(const RegionType &region) override
void SetSpacing(const double *spacing) override
void SetDirection(const DirectionType &direction) override
TVector TransformPhysicalVectorToLocalVector(const TVector &inputGradient) const
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, Self::ImageDimension > &inputGradient, FixedArray< TCoordRep, Self::ImageDimension > &outputGradient) const
ModifiedTimeType GetMTime() const override
void UpdateOutputData() override
Point< TCoordRep, Self::ImageDimension > TransformIndexToPhysicalPoint(const IndexType &index) const
AccessorType & GetPixelAccessor()
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, Self::ImageDimension > &point) const
const AccessorType & GetPixelAccessor() const
TVector TransformLocalVectorToPhysicalVector(const TVector &inputGradient) const
typename TAccessor::ExternalType PixelType
void Allocate(bool initialize=false) override
void SetPixelAccessor(const AccessorType &accessor)
const RegionType & GetLargestPossibleRegion() const override
void SetRequestedRegionToLargestPossibleRegion() override
TAccessor AccessorType
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TCoordRep, Self::ImageDimension > &index, Point< TCoordRep, Self::ImageDimension > &point) const
void SetOrigin(const PointType origin) override
~ImageAdaptor() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
typename TImage::PixelContainer PixelContainer
void SetSpacing(const float *spacing) override
const OffsetValueType * GetOffsetTable() const
PixelType GetPixel(const IndexType &index) const
void SetLargestPossibleRegion(const RegionType &region) override
ContinuousIndex< TIndexRep, TImage::ImageDimension > TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, TImage::ImageDimension > &point) const
const InternalPixelType * GetBufferPointer() const
virtual void SetImage(TImage *)
void UpdateAccessor(T *)
const RegionType & GetRequestedRegion() const override
IndexType ComputeIndex(OffsetValueType offset) const
const SpacingType & GetSpacing() const override
void Graft(const DataObject *data) override
void Update() override
typename TImage::PixelContainerConstPointer PixelContainerConstPointer
InternalPixelType * GetBufferPointer()
virtual void Graft(const Self *imgData)
const RegionType & GetBufferedRegion() const override
IndexType TransformPhysicalPointToIndex(const Point< TCoordRep, Self::ImageDimension > &point) const
typename TAccessor::InternalType InternalPixelType
Base class for templated image classes.
Definition: itkImageBase.h:115
An image region represents a structured region of data.
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for most ITK classes.
Definition: itkObject.h:62
Templated n-dimensional vector image class.
Implements a weak reference to an object.
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
long IndexValueType
Definition: itkIntTypes.h:90
SizeValueType ModifiedTimeType
Definition: itkIntTypes.h:102
unsigned long SizeValueType
Definition: itkIntTypes.h:83
long OffsetValueType
Definition: itkIntTypes.h:94