ITK  6.0.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
228 operator[](const IndexType & index) const
229 {
230 return m_PixelAccessor.Get(m_Image->GetPixel(index));
231 }
232
234 const OffsetValueType *
236
240
243 using PixelContainer = typename TImage::PixelContainer;
244 using PixelContainerPointer = typename TImage::PixelContainerPointer;
245 using PixelContainerConstPointer = typename TImage::PixelContainerConstPointer;
246
250 {
251 return m_Image->GetPixelContainer();
252 }
253
254 const PixelContainer *
256 {
257 return m_Image->GetPixelContainer();
258 }
259
262 void
264
275 virtual void
276 Graft(const Self * imgData);
277
278#ifndef ITK_FUTURE_LEGACY_REMOVE
280 using InternalPixelPointerType [[deprecated("Please just use `InternalPixelType *` instead!")]] = InternalPixelType *;
281#endif
282
287
288 const InternalPixelType *
290
292 void
293 SetSpacing(const SpacingType & spacing) override;
294
295 void
296 SetSpacing(const double * spacing /*[ImageDimension]*/) override;
297
298 void
299 SetSpacing(const float * spacing /*[ImageDimension]*/) override;
300
304 const SpacingType &
305 GetSpacing() const override;
306
310 const PointType &
311 GetOrigin() const override;
312
314 void
315 SetOrigin(const PointType origin) override;
316
317 void
318 SetOrigin(const double * origin /*[ImageDimension]*/) override;
319
320 void
321 SetOrigin(const float * origin /*[ImageDimension]*/) override;
322
324 void
325 SetDirection(const DirectionType & direction) override;
326
330 const DirectionType &
331 GetDirection() const override;
332
334 virtual void
335 SetImage(TImage *);
336
338 void
339 Modified() const override;
340
343 GetMTime() const override;
344
348 {
349 return m_PixelAccessor;
350 }
351
353 const AccessorType &
355 {
356 return m_PixelAccessor;
357 }
358
360 void
362 {
363 m_PixelAccessor = accessor;
364 }
365
367 void
368 Update() override;
369
370 void
371 CopyInformation(const DataObject * data) override;
372
375 void
377
378 void
380
381 void
383
384 void
386
387 bool
389
391 template <typename TIndexRep, typename TCoordinate>
394 {
395 return m_Image->template TransformPhysicalPointToContinuousIndex<TIndexRep>(point);
396 }
397
406 template <typename TCoordinate>
407 ITK_NODISCARD("Call the overload which has the point as the only parameter and returns the index")
408 bool TransformPhysicalPointToContinuousIndex(const Point<TCoordinate, Self::ImageDimension> & point,
409 ContinuousIndex<TCoordinate, Self::ImageDimension> & index) const
410 {
411 return m_Image->TransformPhysicalPointToContinuousIndex(point, index);
412 }
413
415 template <typename TCoordinate>
416 [[nodiscard]] IndexType
418 {
419 return m_Image->TransformPhysicalPointToIndex(point);
420 }
421
430 template <typename TCoordinate>
431 ITK_NODISCARD("Call the overload which has the point as the only parameter and returns the index")
432 bool TransformPhysicalPointToIndex(const Point<TCoordinate, Self::ImageDimension> & point, IndexType & index) const
433 {
434 return m_Image->TransformPhysicalPointToIndex(point, index);
435 }
436
441 template <typename TCoordinate>
442 void
445 {
446 m_Image->TransformContinuousIndexToPhysicalPoint(index, point);
447 }
448
450 template <typename TCoordinate, typename TIndexRep>
453 {
454 return m_Image->template TransformContinuousIndexToPhysicalPoint<TIndexRep>(index);
455 }
456
462 template <typename TCoordinate>
463 void
465 {
466 m_Image->TransformIndexToPhysicalPoint(index, point);
467 }
468
470 template <typename TCoordinate>
473 {
474 return m_Image->template TransformIndexToPhysicalPoint<TCoordinate>(index);
475 }
476
477 template <typename TCoordinate>
478 void
481 {
482 m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
483 }
484
485 template <typename TVector>
486 [[nodiscard]] TVector
487 TransformLocalVectorToPhysicalVector(const TVector & inputGradient) const
488 {
489 TVector outputGradient;
490 TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
491 return outputGradient;
492 }
493
494 template <typename TCoordinate>
495 void
498 {
499 m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
500 }
501
502 template <typename TVector>
503 [[nodiscard]] TVector
504 TransformPhysicalVectorToLocalVector(const TVector & inputGradient) const
505 {
506 TVector outputGradient;
507 TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
508 return outputGradient;
509 }
510
511protected:
513 ~ImageAdaptor() override = default;
514 void
515 PrintSelf(std::ostream & os, Indent indent) const override;
516 void
517 Graft(const DataObject * data) override;
518 using Superclass::Graft;
519
520private:
521 // a specialized method to update PixelAccessors for VectorImages,
522 // to have the correct vector length of the image.
523 template <typename TPixelType>
524 void
526 {
527 this->m_PixelAccessor.SetVectorLength(this->m_Image->GetNumberOfComponentsPerPixel());
528 }
529
530 // The other image types don't expect an accessor which needs any updates
531 template <typename T>
532 void
533 UpdateAccessor(T * itkNotUsed(dummy))
534 {}
535
536 // Adapted image, most of the calls to ImageAdaptor
537 // will be delegated to this image
538 typename TImage::Pointer m_Image{};
539
540 // Data accessor object,
541 // it converts the presentation of a pixel
542 AccessorType m_PixelAccessor{};
543};
544} // end namespace itk
545
546#ifndef ITK_MANUAL_INSTANTIATION
547# include "itkImageAdaptor.hxx"
548#endif
549
550#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
typename SizeType::SizeValueType SizeValueType
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordinate, Self::ImageDimension > &inputGradient, FixedArray< TCoordinate, Self::ImageDimension > &outputGradient) const
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
typename IndexType::IndexValueType IndexValueType
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordinate, Self::ImageDimension > &inputGradient, FixedArray< TCoordinate, Self::ImageDimension > &outputGradient) const
const DirectionType & GetDirection() const override
bool VerifyRequestedRegion() override
void UpdateOutputInformation() override
void SetSpacing(const SpacingType &spacing) override
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordinate, Self::ImageDimension > &point) const
void SetBufferedRegion(const RegionType &region) override
void SetSpacing(const double *spacing) override
Point< TCoordinate, TImage::ImageDimension > TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, Self::ImageDimension > &index) const
void SetDirection(const DirectionType &direction) override
TVector TransformPhysicalVectorToLocalVector(const TVector &inputGradient) const
ModifiedTimeType GetMTime() const override
void UpdateOutputData() override
AccessorType & GetPixelAccessor()
const AccessorType & GetPixelAccessor() const
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TCoordinate, Self::ImageDimension > &index, Point< TCoordinate, Self::ImageDimension > &point) 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
Point< TCoordinate, Self::ImageDimension > TransformIndexToPhysicalPoint(const IndexType &index) const
TAccessor AccessorType
IndexType TransformPhysicalPointToIndex(const Point< TCoordinate, 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
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
ContinuousIndex< TIndexRep, TImage::ImageDimension > TransformPhysicalPointToContinuousIndex(const Point< TCoordinate, TImage::ImageDimension > &point) const
typename TImage::PixelContainerConstPointer PixelContainerConstPointer
InternalPixelType * GetBufferPointer()
virtual void Graft(const Self *imgData)
const RegionType & GetBufferedRegion() const override
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:93
SizeValueType ModifiedTimeType
Definition: itkIntTypes.h:105
unsigned long SizeValueType
Definition: itkIntTypes.h:86
long OffsetValueType
Definition: itkIntTypes.h:97