 |
ITK
5.4.0
Insight Toolkit
|
Go to the documentation of this file.
18 #ifndef itkImageAdaptor_h
19 #define itkImageAdaptor_h
26 template <
typename TPixelType,
unsigned int VImageDimension>
55 template <
typename TImage,
typename TAccessor>
65 static constexpr
unsigned int ImageDimension = TImage::ImageDimension;
99 using AccessorFunctorType =
typename InternalImageType::AccessorFunctorType::template Rebind<Self>::Type;
110 using typename Superclass::OffsetType;
119 using typename Superclass::SpacingType;
137 template <
typename UPixelType,
unsigned int UImageDimension = TImage::ImageDimension>
143 template <
typename UPixelType,
unsigned int VUImageDimension = TImage::ImageDimension>
154 SetLargestPossibleRegion(
const RegionType & region)
override;
160 SetBufferedRegion(
const RegionType & region)
override;
166 SetRequestedRegion(
const RegionType & region)
override;
173 SetRequestedRegion(
const DataObject * data)
override;
182 GetRequestedRegion()
const override;
193 GetLargestPossibleRegion()
const override;
201 GetBufferedRegion()
const override;
205 Allocate(
bool initialize =
false)
override;
210 Initialize()
override;
216 m_PixelAccessor.Set(m_Image->GetPixel(index), value);
223 return m_PixelAccessor.Get(m_Image->GetPixel(index));
231 GetOffsetTable()
const;
247 return m_Image->GetPixelContainer();
250 const PixelContainer *
253 return m_Image->GetPixelContainer();
259 SetPixelContainer(PixelContainer * container);
272 Graft(
const Self * imgData);
283 GetBufferPointer()
const;
290 SetSpacing(
const double * spacing )
override;
293 SetSpacing(
const float * spacing )
override;
299 GetSpacing()
const override;
305 GetOrigin()
const override;
309 SetOrigin(
const PointType origin)
override;
312 SetOrigin(
const double * origin )
override;
315 SetOrigin(
const float * origin )
override;
325 GetDirection()
const override;
333 Modified()
const override;
337 GetMTime()
const override;
343 return m_PixelAccessor;
350 return m_PixelAccessor;
357 m_PixelAccessor = accessor;
365 CopyInformation(
const DataObject * data)
override;
370 UpdateOutputInformation()
override;
373 SetRequestedRegionToLargestPossibleRegion()
override;
376 PropagateRequestedRegion()
override;
379 UpdateOutputData()
override;
382 VerifyRequestedRegion()
override;
385 template <
typename TIndexRep,
typename TCoordRep>
389 return m_Image->template TransformPhysicalPointToContinuousIndex<TIndexRep>(
point);
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,
405 return m_Image->TransformPhysicalPointToContinuousIndex(
point, index);
409 template <
typename TCoordRep>
413 return m_Image->TransformPhysicalPointToIndex(
point);
424 template <
typename TCoordRep>
425 ITK_NODISCARD(
"Call the overload which has the point as the only parameter and returns the index")
428 return m_Image->TransformPhysicalPointToIndex(
point, index);
435 template <
typename TCoordRep>
440 m_Image->TransformContinuousIndexToPhysicalPoint(index,
point);
444 template <
typename TCoordRep,
typename TIndexRep>
448 return m_Image->template TransformContinuousIndexToPhysicalPoint<TIndexRep>(index);
456 template <
typename TCoordRep>
460 m_Image->TransformIndexToPhysicalPoint(index,
point);
464 template <
typename TCoordRep>
468 return m_Image->template TransformIndexToPhysicalPoint<TCoordRep>(index);
471 template <
typename TCoordRep>
476 m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
479 template <
typename TVector>
480 [[nodiscard]] TVector
483 TVector outputGradient;
484 TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
485 return outputGradient;
488 template <
typename TCoordRep>
493 m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
496 template <
typename TVector>
497 [[nodiscard]] TVector
500 TVector outputGradient;
501 TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
502 return outputGradient;
509 PrintSelf(std::ostream & os,
Indent indent)
const override;
512 using Superclass::Graft;
517 template <
typename TPixelType>
521 this->m_PixelAccessor.SetVectorLength(this->m_Image->GetNumberOfComponentsPerPixel());
525 template <
typename T>
540 #ifndef ITK_MANUAL_INSTANTIATION
541 # include "itkImageAdaptor.hxx"
PixelType GetPixel(const IndexType &index) const
typename OffsetType::OffsetValueType OffsetValueType
typename TImage::PixelContainerConstPointer PixelContainerConstPointer
SmartPointer< Self > Pointer
ImageBaseType::DirectionType DirectionType
SizeValueType ModifiedTimeType
typename TImage::PixelContainer PixelContainer
Base class for templated image classes.
An image region represents a structured region of data.
TVector TransformLocalVectorToPhysicalVector(const TVector &inputGradient) const
typename SizeType::SizeValueType SizeValueType
typename Accessor::AddPixelAccessor< TImage::PixelType > ::InternalType InternalPixelType
Templated n-dimensional vector image class.
ImageBaseType::PointType PointType
ImageBaseType::SizeType SizeType
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, Self::ImageDimension > &inputGradient, FixedArray< TCoordRep, Self::ImageDimension > &outputGradient) const
InternalPixelType * InternalPixelPointerType
Control indentation during Print() invocation.
const PixelContainer * GetPixelContainer() const
typename IndexType::IndexValueType IndexValueType
PixelType operator[](const IndexType &index) const
Give access to partial aspects of voxels from an Image.
ImageBaseType::IndexType IndexType
ContinuousIndex< TIndexRep, TImage::ImageDimension > TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, TImage::ImageDimension > &point) const
IndexType TransformPhysicalPointToIndex(const Point< TCoordRep, Self::ImageDimension > &point) const
void UpdateAccessor(typename itk::VectorImage< TPixelType, ImageDimension > *)
TVector TransformPhysicalVectorToLocalVector(const TVector &inputGradient) const
*par Constraints *The filter requires an 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
ImageBaseType::RegionType RegionType
typename InternalImageType::AccessorFunctorType::template Rebind< Self >::Type AccessorFunctorType
void SetPixelAccessor(const AccessorType &accessor)
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, Self::ImageDimension > &point) const
Simulate a standard C array with copy semantics.
typename TImage::PixelContainerPointer PixelContainerPointer
typename Accessor::AddPixelAccessor< TImage::PixelType > ::ExternalType PixelType
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, Self::ImageDimension > &inputGradient, FixedArray< TCoordRep, Self::ImageDimension > &outputGradient) const
Implements a weak reference to an object.
Point< TCoordRep, Self::ImageDimension > TransformIndexToPhysicalPoint(const IndexType &index) const
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TCoordRep, Self::ImageDimension > &index, Point< TCoordRep, Self::ImageDimension > &point) const
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
A templated class holding a point in n-Dimensional image space.
Base class for most ITK classes.
void SetPixel(const IndexType &index, const PixelType &value)
Templated n-dimensional image class.
AccessorType & GetPixelAccessor()
const AccessorType & GetPixelAccessor() const
PixelContainerPointer GetPixelContainer()
unsigned long SizeValueType
Base class for all data objects in ITK.
Point< TCoordRep, TImage::ImageDimension > TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, Self::ImageDimension > &index) const