ITK  5.2.0
Insight Toolkit
Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
itk::VectorImage< TPixel, VImageDimension > Class Template Reference

#include <itkImageAlgorithm.h>

+ Inheritance diagram for itk::VectorImage< TPixel, VImageDimension >:
+ Collaboration diagram for itk::VectorImage< TPixel, VImageDimension >:

Classes

struct  Rebind
 

Public Types

using AccessorFunctorType = DefaultVectorPixelAccessorFunctor< Self >
 
using AccessorType = DefaultVectorPixelAccessor< InternalPixelType >
 
using ConstPointer = SmartPointer< const Self >
 
using ConstWeakPointer = WeakPointer< const Self >
 
using DirectionType = typename Superclass::DirectionType
 
using IndexType = typename Superclass::IndexType
 
using IndexValueType = typename Superclass::IndexValueType
 
using InternalPixelType = TPixel
 
using IOPixelType = InternalPixelType
 
using NeighborhoodAccessorFunctorType = VectorImageNeighborhoodAccessorFunctor< Self >
 
using OffsetType = typename Superclass::OffsetType
 
using OffsetValueType = typename Superclass::OffsetValueType
 
using PixelContainer = ImportImageContainer< SizeValueType, InternalPixelType >
 
using PixelContainerConstPointer = typename PixelContainer::ConstPointer
 
using PixelContainerPointer = typename PixelContainer::Pointer
 
using PixelType = VariableLengthVector< TPixel >
 
using Pointer = SmartPointer< Self >
 
using PointType = typename Superclass::PointType
 
template<typename UPixelType , unsigned int NUImageDimension = VImageDimension>
using RebindImageType = typename Rebind< UPixelType, NUImageDimension >::Type
 
using RegionType = typename Superclass::RegionType
 
using Self = VectorImage
 
using SizeType = typename Superclass::SizeType
 
using SpacingType = typename Superclass::SpacingType
 
using Superclass = ImageBase< VImageDimension >
 
using ValueType = PixelType
 
using VectorLengthType = unsigned int
 
- Public Types inherited from itk::ImageBase< VImageDimension >
using ConstPointer = SmartPointer< const Self >
 
using DirectionType = Matrix< SpacePrecisionType, VImageDimension, VImageDimension >
 
using ImageDimensionType = unsigned int
 
using IndexType = Index< VImageDimension >
 
using IndexValueType = typename IndexType::IndexValueType
 
using OffsetType = Offset< VImageDimension >
 
using OffsetValueType = typename OffsetType::OffsetValueType
 
using Pointer = SmartPointer< Self >
 
using PointType = Point< PointValueType, VImageDimension >
 
using PointValueType = SpacePrecisionType
 
using RegionType = ImageRegion< VImageDimension >
 
using Self = ImageBase
 
using SizeType = Size< VImageDimension >
 
using SizeValueType = typename SizeType::SizeValueType
 
using SpacingType = Vector< SpacingValueType, VImageDimension >
 
using SpacingValueType = SpacePrecisionType
 
using Superclass = DataObject
 
- Public Types inherited from itk::DataObject
using ConstPointer = SmartPointer< const Self >
 
using DataObjectIdentifierType = std::string
 
using DataObjectPointerArraySizeType = std::vector< Pointer >::size_type
 
using Pointer = SmartPointer< Self >
 
using Self = DataObject
 
using Superclass = Object
 
- Public Types inherited from itk::Object
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = Object
 
using Superclass = LightObject
 
- Public Types inherited from itk::LightObject
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = LightObject
 

Public Member Functions

void Allocate (bool UseDefaultConstructor=false) override
 
virtual ::itk::LightObject::Pointer CreateAnother () const
 
void FillBuffer (const PixelType &value)
 
virtual const char * GetNameOfClass () const
 
PixelType GetPixel (const IndexType &index)
 
const PixelType GetPixel (const IndexType &index) const
 
void Initialize () override
 
PixelType operator[] (const IndexType &index)
 
const PixelType operator[] (const IndexType &index) const
 
void SetPixel (const IndexType &index, const PixelType &value)
 
- Public Member Functions inherited from itk::ImageBase< VImageDimension >
virtual void SetOrigin (PointType _arg)
 
virtual void SetOrigin (const double origin[VImageDimension])
 
virtual void SetOrigin (const float origin[VImageDimension])
 
virtual void SetDirection (const DirectionType &direction)
 
virtual const DirectionTypeGetDirection () const
 
virtual const DirectionTypeGetInverseDirection () const
 
virtual const SpacingTypeGetSpacing () const
 
virtual const PointTypeGetOrigin () const
 
virtual void SetLargestPossibleRegion (const RegionType &region)
 
virtual const RegionTypeGetLargestPossibleRegion () const
 
virtual void SetBufferedRegion (const RegionType &region)
 
virtual const RegionTypeGetBufferedRegion () const
 
virtual void SetRequestedRegion (const RegionType &region)
 
void SetRequestedRegion (const DataObject *data) override
 
virtual const RegionTypeGetRequestedRegion () const
 
virtual void SetRegions (const RegionType &region)
 
virtual void SetRegions (const SizeType &size)
 
const OffsetValueTypeGetOffsetTable () const
 
OffsetValueType ComputeOffset (const IndexType &ind) const
 
IndexType ComputeIndex (OffsetValueType offset) const
 
virtual void SetSpacing (const SpacingType &spacing)
 
virtual void SetSpacing (const double spacing[VImageDimension])
 
virtual void SetSpacing (const float spacing[VImageDimension])
 
template<typename TCoordRep >
IndexType TransformPhysicalPointToIndex (const Point< TCoordRep, VImageDimension > &point) const
 
template<typename TCoordRep >
bool TransformPhysicalPointToIndex (const Point< TCoordRep, VImageDimension > &point, IndexType &index) const
 
template<typename TIndexRep , typename TCoordRep >
ContinuousIndex< TIndexRep, VImageDimension > TransformPhysicalPointToContinuousIndex (const Point< TCoordRep, VImageDimension > &point) const
 
template<typename TCoordRep , typename TIndexRep >
bool TransformPhysicalPointToContinuousIndex (const Point< TCoordRep, VImageDimension > &point, ContinuousIndex< TIndexRep, VImageDimension > &index) const
 
template<typename TCoordRep , typename TIndexRep >
void TransformContinuousIndexToPhysicalPoint (const ContinuousIndex< TIndexRep, VImageDimension > &index, Point< TCoordRep, VImageDimension > &point) const
 
template<typename TCoordRep , typename TIndexRep >
Point< TCoordRep, VImageDimension > TransformContinuousIndexToPhysicalPoint (const ContinuousIndex< TIndexRep, VImageDimension > &index) const
 
template<typename TCoordRep >
void TransformIndexToPhysicalPoint (const IndexType &index, Point< TCoordRep, VImageDimension > &point) const
 
template<typename TCoordRep >
Point< TCoordRep, VImageDimension > TransformIndexToPhysicalPoint (const IndexType &index) const
 
template<typename TCoordRep >
void TransformLocalVectorToPhysicalVector (const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
 
template<typename TVector >
TVector TransformLocalVectorToPhysicalVector (const TVector &inputGradient) const
 
template<typename TCoordRep >
void TransformPhysicalVectorToLocalVector (const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
 
template<typename TVector >
TVector TransformPhysicalVectorToLocalVector (const TVector &inputGradient) const
 
void CopyInformation (const DataObject *data) override
 
void UpdateOutputInformation () override
 
void UpdateOutputData () override
 
void SetRequestedRegionToLargestPossibleRegion () override
 
bool RequestedRegionIsOutsideOfTheBufferedRegion () override
 
bool VerifyRequestedRegion () override
 
- Public Member Functions inherited from itk::DataObject
void DisconnectPipeline ()
 
bool GetDataReleased () const
 
virtual const bool & GetReleaseDataFlag () const
 
SmartPointer< ProcessObjectGetSource () const
 
DataObjectPointerArraySizeType GetSourceOutputIndex () const
 
const DataObjectIdentifierTypeGetSourceOutputName () const
 
virtual void PropagateRequestedRegion ()
 
void ReleaseData ()
 
virtual void ReleaseDataFlagOff ()
 
virtual void ReleaseDataFlagOn ()
 
virtual void ResetPipeline ()
 
void SetReleaseDataFlag (bool flag)
 
bool ShouldIReleaseData () const
 
virtual void Update ()
 
void SetPipelineMTime (ModifiedTimeType time)
 
virtual const ModifiedTimeTypeGetPipelineMTime () const
 
virtual ModifiedTimeType GetUpdateMTime () const
 
virtual void SetRealTimeStamp (RealTimeStamp _arg)
 
virtual const RealTimeStampGetRealTimeStamp () const
 
virtual void PrepareForNewData ()
 
virtual void DataHasBeenGenerated ()
 
- Public Member Functions inherited from itk::Object
unsigned long AddObserver (const EventObject &event, Command *)
 
unsigned long AddObserver (const EventObject &event, Command *) const
 
unsigned long AddObserver (const EventObject &event, std::function< void(const EventObject &)> function) const
 
virtual void DebugOff () const
 
virtual void DebugOn () const
 
CommandGetCommand (unsigned long tag)
 
bool GetDebug () const
 
MetaDataDictionaryGetMetaDataDictionary ()
 
const MetaDataDictionaryGetMetaDataDictionary () const
 
virtual ModifiedTimeType GetMTime () const
 
virtual const TimeStampGetTimeStamp () const
 
bool HasObserver (const EventObject &event) const
 
void InvokeEvent (const EventObject &)
 
void InvokeEvent (const EventObject &) const
 
virtual void Modified () const
 
void Register () const override
 
void RemoveAllObservers ()
 
void RemoveObserver (unsigned long tag)
 
void SetDebug (bool debugFlag) const
 
void SetReferenceCount (int) override
 
void UnRegister () const noexcept override
 
void SetMetaDataDictionary (const MetaDataDictionary &rhs)
 
void SetMetaDataDictionary (MetaDataDictionary &&rrhs)
 
virtual void SetObjectName (std::string _arg)
 
virtual const std::string & GetObjectName () const
 
- Public Member Functions inherited from itk::LightObject
Pointer Clone () const
 
virtual void Delete ()
 
virtual int GetReferenceCount () const
 
void Print (std::ostream &os, Indent indent=0) const
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from itk::ImageBase< VImageDimension >
static unsigned int GetImageDimension ()
 
static Pointer New ()
 
- Static Public Member Functions inherited from itk::DataObject
static bool GetGlobalReleaseDataFlag ()
 
static void GlobalReleaseDataFlagOff ()
 
static void GlobalReleaseDataFlagOn ()
 
static Pointer New ()
 
static void SetGlobalReleaseDataFlag (bool val)
 
- Static Public Member Functions inherited from itk::Object
static bool GetGlobalWarningDisplay ()
 
static void GlobalWarningDisplayOff ()
 
static void GlobalWarningDisplayOn ()
 
static Pointer New ()
 
static void SetGlobalWarningDisplay (bool val)
 
- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 

Static Public Attributes

static constexpr unsigned int ImageDimension = VImageDimension
 
- Static Public Attributes inherited from itk::ImageBase< VImageDimension >
static constexpr ImageDimensionType ImageDimension = VImageDimension
 
VectorLengthType m_VectorLength { 0 }
 
PixelContainerPointer m_Buffer
 
InternalPixelTypeGetBufferPointer ()
 
const InternalPixelTypeGetBufferPointer () const
 
PixelContainerGetPixelContainer ()
 
const PixelContainerGetPixelContainer () const
 
void SetPixelContainer (PixelContainer *container)
 
virtual void Graft (const Self *image)
 
AccessorType GetPixelAccessor ()
 
const AccessorType GetPixelAccessor () const
 
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor ()
 
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor () const
 
virtual void SetVectorLength (VectorLengthType _arg)
 
virtual const VectorLengthTypeGetVectorLength () const
 
unsigned int GetNumberOfComponentsPerPixel () const override
 
void SetNumberOfComponentsPerPixel (unsigned int n) override
 
 VectorImage ()
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~VectorImage () override=default
 
void Graft (const DataObject *data) override
 

Additional Inherited Members

- Protected Member Functions inherited from itk::ImageBase< VImageDimension >
 ImageBase ()
 
 ~ImageBase () override=default
 
void ComputeOffsetTable ()
 
virtual void ComputeIndexToPhysicalPointMatrices ()
 
virtual void InitializeBufferedRegion ()
 
OffsetValueType FastComputeOffset (const IndexType &ind) const
 
IndexType FastComputeIndex (OffsetValueType offset) const
 
- Protected Member Functions inherited from itk::DataObject
 DataObject ()
 
 ~DataObject () override
 
virtual void PropagateResetPipeline ()
 
- Protected Member Functions inherited from itk::Object
 Object ()
 
 ~Object () override
 
bool PrintObservers (std::ostream &os, Indent indent) const
 
virtual void SetTimeStamp (const TimeStamp &timeStamp)
 
- Protected Member Functions inherited from itk::LightObject
virtual LightObject::Pointer InternalClone () const
 
 LightObject ()
 
virtual void PrintHeader (std::ostream &os, Indent indent) const
 
virtual void PrintTrailer (std::ostream &os, Indent indent) const
 
virtual ~LightObject ()
 
- Protected Attributes inherited from itk::ImageBase< VImageDimension >
SpacingType m_Spacing
 
PointType m_Origin
 
DirectionType m_Direction
 
DirectionType m_InverseDirection
 
DirectionType m_IndexToPhysicalPoint
 
DirectionType m_PhysicalPointToIndex
 
- Protected Attributes inherited from itk::LightObject
std::atomic< int > m_ReferenceCount
 

Detailed Description

template<typename TPixel, unsigned int VImageDimension = 3>
class itk::VectorImage< TPixel, VImageDimension >

Templated n-dimensional vector image class.

This class differs from Image in that it is intended to represent multiple images. Each pixel represents k measurements, each of datatype TPixel. The memory organization of the resulting image is as follows: ... Pi0 Pi1 Pi2 Pi3 P(i+1)0 P(i+1)1 P(i+1)2 P(i+1)3 P(i+2)0 ... where Pi0 represents the 0th measurement of the pixel at index i.

Conceptually, a VectorImage< TPixel, 3 > is the same as a Image< VariableLengthVector< TPixel >, 3 >. The difference lies in the memory organization. The latter results in a fragmented organization with each location in the Image holding a pointer to an VariableLengthVector holding the actual pixel. The former stores the k pixels instead of a pointer reference, which apart from avoiding fragmentation of memory also avoids storing a 8 bytes of pointer reference for each pixel. The parameter k can be set using SetVectorLength.

The API of the class is such that it returns a pixeltype VariableLengthVector< TPixel > when queried, with the data internally pointing to the buffer. (the container does not manage the memory). Similarly SetPixel calls can be made with VariableLengthVector< TPixel >.

The API of this class is similar to Image.

Caveats:
When using Iterators on this image, you cannot use the it.Value(). You must use Set/Get() methods instead.
Note
This work is part of the National Alliance for Medical Image Computing (NAMIC), funded by the National Institutes of Health through the NIH Roadmap for Medical Research, Grant U54 EB005149.
See also
DefaultVectorPixelAccessor
DefaultVectorPixelAccessorFunctor
VectorImageToImagePixelAccessor
VectorImageToImageAdaptor
Image
ImportImageContainer
ITK Sphinx Examples:
Examples
Examples/Filtering/DiffusionTensor3DReconstructionImageFilter.cxx, Examples/Statistics/BayesianClassifier.cxx, SphinxExamples/src/Core/Common/ApplyCustomOperationToEachPixelInImage/Code.cxx, SphinxExamples/src/Core/Common/CastVectorImageToAnotherType/Code.cxx, SphinxExamples/src/Core/Common/CreateVectorImage/Code.cxx, SphinxExamples/src/Core/Common/NeighborhoodIteratorOnVectorImage/Code.cxx, SphinxExamples/src/Core/ImageAdaptors/ViewComponentVectorImageAsScaleImage/Code.cxx, SphinxExamples/src/Filtering/ImageCompose/CreateVectorImageFromScalarImages/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/ComputerMagInVectorImageToMakeMagImage/Code.cxx, and SphinxExamples/src/Filtering/ImageIntensity/ExtractComponentOfVectorImage/Code.cxx.

Definition at line 29 of file itkImageAlgorithm.h.

Member Typedef Documentation

◆ AccessorFunctorType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::AccessorFunctorType = DefaultVectorPixelAccessorFunctor<Self>

Functor to provide a common API between DefaultPixelAccessor and DefaultVectorPixelAccessor

Definition at line 121 of file itkVectorImage.h.

◆ AccessorType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::AccessorType = DefaultVectorPixelAccessor<InternalPixelType>

Accessor type that convert data between internal and external representations.

Definition at line 117 of file itkVectorImage.h.

◆ ConstPointer

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::ConstPointer = SmartPointer<const Self>

Definition at line 90 of file itkVectorImage.h.

◆ ConstWeakPointer

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::ConstWeakPointer = WeakPointer<const Self>

Definition at line 91 of file itkVectorImage.h.

◆ DirectionType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::DirectionType = typename Superclass::DirectionType

Direction type alias support A matrix of direction cosines.

Definition at line 147 of file itkVectorImage.h.

◆ IndexType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::IndexType = typename Superclass::IndexType

Index type alias support An index is used to access pixel values.

Definition at line 134 of file itkVectorImage.h.

◆ IndexValueType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::IndexValueType = typename Superclass::IndexValueType

Definition at line 135 of file itkVectorImage.h.

◆ InternalPixelType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::InternalPixelType = TPixel

This is the actual pixel type contained in the buffer. Each vector pixel is composed of 'm_VectorLength' contiguous InternalPixelType.

Definition at line 108 of file itkVectorImage.h.

◆ IOPixelType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::IOPixelType = InternalPixelType

Definition at line 113 of file itkVectorImage.h.

◆ NeighborhoodAccessorFunctorType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::NeighborhoodAccessorFunctorType = VectorImageNeighborhoodAccessorFunctor<Self>

Typedef for the functor used to access a neighborhood of pixel pointers.

Definition at line 125 of file itkVectorImage.h.

◆ OffsetType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::OffsetType = typename Superclass::OffsetType

Offset type alias support An offset is used to access pixel values.

Definition at line 138 of file itkVectorImage.h.

◆ OffsetValueType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::OffsetValueType = typename Superclass::OffsetValueType

Offset type alias (relative position between indices)

Definition at line 166 of file itkVectorImage.h.

◆ PixelContainer

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::PixelContainer = ImportImageContainer<SizeValueType, InternalPixelType>

Container used to store pixels in the image.

Definition at line 144 of file itkVectorImage.h.

◆ PixelContainerConstPointer

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::PixelContainerConstPointer = typename PixelContainer::ConstPointer

Definition at line 163 of file itkVectorImage.h.

◆ PixelContainerPointer

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::PixelContainerPointer = typename PixelContainer::Pointer

A pointer to the pixel container.

Definition at line 162 of file itkVectorImage.h.

◆ PixelType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::PixelType = VariableLengthVector<TPixel>

Pixel type alias support Used to declare pixel type in filters or other operations. This is not the actual pixel type contained in the buffer, ie m_Buffer. The image exhibits an external API of an VariableLengthVector< T > and internally stores its data as type T.

Definition at line 103 of file itkVectorImage.h.

◆ Pointer

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::Pointer = SmartPointer<Self>

Definition at line 89 of file itkVectorImage.h.

◆ PointType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::PointType = typename Superclass::PointType

Origin type alias support The origin is the geometric coordinates of the index (0,0).

Definition at line 159 of file itkVectorImage.h.

◆ RebindImageType

template<typename TPixel, unsigned int VImageDimension = 3>
template<typename UPixelType , unsigned int NUImageDimension = VImageDimension>
using itk::VectorImage< TPixel, VImageDimension >::RebindImageType = typename Rebind<UPixelType, NUImageDimension>::Type

Definition at line 201 of file itkVectorImage.h.

◆ RegionType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::RegionType = typename Superclass::RegionType

Region type alias support A region is used to specify a subset of an image.

Definition at line 151 of file itkVectorImage.h.

◆ Self

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::Self = VectorImage

Standard class type aliases

Definition at line 87 of file itkVectorImage.h.

◆ SizeType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::SizeType = typename Superclass::SizeType

Size type alias support A size is used to define region bounds.

Definition at line 141 of file itkVectorImage.h.

◆ SpacingType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::SpacingType = typename Superclass::SpacingType

Spacing type alias support Spacing holds the size of a pixel. The spacing is the geometric distance between image samples.

Definition at line 155 of file itkVectorImage.h.

◆ Superclass

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::Superclass = ImageBase<VImageDimension>

Definition at line 88 of file itkVectorImage.h.

◆ ValueType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::ValueType = PixelType

Typedef alias for PixelType

Definition at line 111 of file itkVectorImage.h.

◆ VectorLengthType

template<typename TPixel, unsigned int VImageDimension = 3>
using itk::VectorImage< TPixel, VImageDimension >::VectorLengthType = unsigned int

Definition at line 168 of file itkVectorImage.h.

Constructor & Destructor Documentation

◆ VectorImage()

template<typename TPixel, unsigned int VImageDimension = 3>
itk::VectorImage< TPixel, VImageDimension >::VectorImage ( )
protected

Length of the "vector pixel"

◆ ~VectorImage()

template<typename TPixel, unsigned int VImageDimension = 3>
itk::VectorImage< TPixel, VImageDimension >::~VectorImage ( )
overrideprotecteddefault

Length of the "vector pixel"

Member Function Documentation

◆ Allocate()

template<typename TPixel, unsigned int VImageDimension = 3>
void itk::VectorImage< TPixel, VImageDimension >::Allocate ( bool  UseDefaultConstructor = false)
overridevirtual

Allocate the image memory. The size of the image must already be set, e.g. by calling SetRegions().

Reimplemented from itk::ImageBase< VImageDimension >.

◆ CreateAnother()

template<typename TPixel, unsigned int VImageDimension = 3>
virtual::itk::LightObject::Pointer itk::VectorImage< TPixel, VImageDimension >::CreateAnother ( ) const
virtual

Create an object from an instance, potentially deferring to a factory. This method allows you to create an instance of an object that is exactly the same type as the referring object. This is useful in cases where an object has been cast back to a base class.

Reimplemented from itk::ImageBase< VImageDimension >.

◆ FillBuffer()

template<typename TPixel, unsigned int VImageDimension = 3>
void itk::VectorImage< TPixel, VImageDimension >::FillBuffer ( const PixelType value)

Fill the image buffer with a value. Be sure to call Allocate() first.

◆ GetBufferPointer() [1/2]

template<typename TPixel, unsigned int VImageDimension = 3>
InternalPixelType* itk::VectorImage< TPixel, VImageDimension >::GetBufferPointer ( )
inline

Return a pointer to the beginning of the buffer. This is used by the image iterator class.

Definition at line 288 of file itkVectorImage.h.

◆ GetBufferPointer() [2/2]

template<typename TPixel, unsigned int VImageDimension = 3>
const InternalPixelType* itk::VectorImage< TPixel, VImageDimension >::GetBufferPointer ( ) const
inline

Length of the "vector pixel"

Definition at line 293 of file itkVectorImage.h.

◆ GetNameOfClass()

template<typename TPixel, unsigned int VImageDimension = 3>
virtual const char* itk::VectorImage< TPixel, VImageDimension >::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from itk::ImageBase< VImageDimension >.

◆ GetNeighborhoodAccessor() [1/2]

template<typename TPixel, unsigned int VImageDimension = 3>
NeighborhoodAccessorFunctorType itk::VectorImage< TPixel, VImageDimension >::GetNeighborhoodAccessor ( )
inline

Return the NeighborhoodAccessor functor

Definition at line 347 of file itkVectorImage.h.

◆ GetNeighborhoodAccessor() [2/2]

template<typename TPixel, unsigned int VImageDimension = 3>
const NeighborhoodAccessorFunctorType itk::VectorImage< TPixel, VImageDimension >::GetNeighborhoodAccessor ( ) const
inline

Return the NeighborhoodAccessor functor

Definition at line 354 of file itkVectorImage.h.

◆ GetNumberOfComponentsPerPixel()

template<typename TPixel, unsigned int VImageDimension = 3>
unsigned int itk::VectorImage< TPixel, VImageDimension >::GetNumberOfComponentsPerPixel ( ) const
overridevirtual

Get/Set the number of components each pixel has, ie the VectorLength

Reimplemented from itk::ImageBase< VImageDimension >.

◆ GetPixel() [1/2]

template<typename TPixel, unsigned int VImageDimension = 3>
PixelType itk::VectorImage< TPixel, VImageDimension >::GetPixel ( const IndexType index)
inline

Get a "reference" to a pixel. This result cannot be used as an lvalue because the pixel is converted on the fly to a VariableLengthVector.

To use the results to modify this image, return value optimization must be relied upon.

For efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 259 of file itkVectorImage.h.

◆ GetPixel() [2/2]

template<typename TPixel, unsigned int VImageDimension = 3>
const PixelType itk::VectorImage< TPixel, VImageDimension >::GetPixel ( const IndexType index) const
inline

Get a pixel (read only version).

For efficiency, this function does not check that the image has actually been allocated yet. Note that the method returns a pixel on the stack.

Definition at line 240 of file itkVectorImage.h.

◆ GetPixelAccessor() [1/2]

template<typename TPixel, unsigned int VImageDimension = 3>
AccessorType itk::VectorImage< TPixel, VImageDimension >::GetPixelAccessor ( )
inline

Return the Pixel Accessor object

Definition at line 333 of file itkVectorImage.h.

◆ GetPixelAccessor() [2/2]

template<typename TPixel, unsigned int VImageDimension = 3>
const AccessorType itk::VectorImage< TPixel, VImageDimension >::GetPixelAccessor ( ) const
inline

Return the Pixel Accesor object

Definition at line 340 of file itkVectorImage.h.

◆ GetPixelContainer() [1/2]

template<typename TPixel, unsigned int VImageDimension = 3>
PixelContainer* itk::VectorImage< TPixel, VImageDimension >::GetPixelContainer ( )
inline

Return a pointer to the container.

Definition at line 301 of file itkVectorImage.h.

◆ GetPixelContainer() [2/2]

template<typename TPixel, unsigned int VImageDimension = 3>
const PixelContainer* itk::VectorImage< TPixel, VImageDimension >::GetPixelContainer ( ) const
inline

Return a pointer to the container.

Definition at line 308 of file itkVectorImage.h.

◆ GetVectorLength()

template<typename TPixel, unsigned int VImageDimension = 3>
virtual const VectorLengthType& itk::VectorImage< TPixel, VImageDimension >::GetVectorLength ( ) const
virtual

Length of the "vector pixel"

◆ Graft() [1/2]

template<typename TPixel, unsigned int VImageDimension = 3>
void itk::VectorImage< TPixel, VImageDimension >::Graft ( const DataObject data)
overrideprotectedvirtual

Length of the "vector pixel"

Reimplemented from itk::ImageBase< VImageDimension >.

◆ Graft() [2/2]

template<typename TPixel, unsigned int VImageDimension = 3>
virtual void itk::VectorImage< TPixel, VImageDimension >::Graft ( const Self image)
virtual

Graft the data and information from one image to another. This is a convenience method to setup a second image with all the meta information of another image and use the same pixel container. Note that this method is different than just using two SmartPointers to the same image since separate DataObjects are still maintained. This method is similar to ImageSource::GraftOutput(). The implementation in ImageBase simply calls CopyInformation() and copies the region ivars. The implementation here refers to the superclass' implementation and then copies over the pixel container.

Reimplemented from itk::ImageBase< VImageDimension >.

◆ Initialize()

template<typename TPixel, unsigned int VImageDimension = 3>
void itk::VectorImage< TPixel, VImageDimension >::Initialize ( )
overridevirtual

Restore the data object to its initial state. This means releasing memory.

Reimplemented from itk::ImageBase< VImageDimension >.

◆ New()

template<typename TPixel, unsigned int VImageDimension = 3>
static Pointer itk::VectorImage< TPixel, VImageDimension >::New ( )
static

Method for creation through the object factory.

◆ operator[]() [1/2]

template<typename TPixel, unsigned int VImageDimension = 3>
PixelType itk::VectorImage< TPixel, VImageDimension >::operator[] ( const IndexType index)
inline

Access a pixel. This result cannot be used as an lvalue because the pixel is converted on the fly to a VariableLengthVector.

To use the results to modify this image, return value optimization must be relied upon.

For efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 277 of file itkVectorImage.h.

◆ operator[]() [2/2]

template<typename TPixel, unsigned int VImageDimension = 3>
const PixelType itk::VectorImage< TPixel, VImageDimension >::operator[] ( const IndexType index) const
inline

Access a pixel.

For efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 283 of file itkVectorImage.h.

◆ PrintSelf()

template<typename TPixel, unsigned int VImageDimension = 3>
void itk::VectorImage< TPixel, VImageDimension >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
overrideprotectedvirtual

Length of the "vector pixel"

Reimplemented from itk::ImageBase< VImageDimension >.

◆ SetNumberOfComponentsPerPixel()

template<typename TPixel, unsigned int VImageDimension = 3>
void itk::VectorImage< TPixel, VImageDimension >::SetNumberOfComponentsPerPixel ( unsigned int  n)
overridevirtual

Length of the "vector pixel"

Reimplemented from itk::ImageBase< VImageDimension >.

◆ SetPixel()

template<typename TPixel, unsigned int VImageDimension = 3>
void itk::VectorImage< TPixel, VImageDimension >::SetPixel ( const IndexType index,
const PixelType value 
)
inline

Set a pixel value.

Allocate() needs to have been called first – for efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 224 of file itkVectorImage.h.

◆ SetPixelContainer()

template<typename TPixel, unsigned int VImageDimension = 3>
void itk::VectorImage< TPixel, VImageDimension >::SetPixelContainer ( PixelContainer container)

Set the container to use. Note that this does not cause the DataObject to be modified.

◆ SetVectorLength()

template<typename TPixel, unsigned int VImageDimension = 3>
virtual void itk::VectorImage< TPixel, VImageDimension >::SetVectorLength ( VectorLengthType  _arg)
virtual

Set/Get macros for the length of each vector in the vector image

Member Data Documentation

◆ ImageDimension

template<typename TPixel, unsigned int VImageDimension = 3>
constexpr unsigned int itk::VectorImage< TPixel, VImageDimension >::ImageDimension = VImageDimension
staticconstexpr

Dimension of the image. This constant is used by functions that are templated over image type (as opposed to being templated over pixel type and dimension) when they need compile time access to the dimension of the image.

Definition at line 131 of file itkVectorImage.h.

◆ m_Buffer

template<typename TPixel, unsigned int VImageDimension = 3>
PixelContainerPointer itk::VectorImage< TPixel, VImageDimension >::m_Buffer
private

Memory for the current buffer.

Definition at line 386 of file itkVectorImage.h.

◆ m_VectorLength

template<typename TPixel, unsigned int VImageDimension = 3>
VectorLengthType itk::VectorImage< TPixel, VImageDimension >::m_VectorLength { 0 }
private

Length of the "vector pixel"

Definition at line 383 of file itkVectorImage.h.


The documentation for this class was generated from the following files: