ITK
4.13.0
Insight Segmentation and Registration Toolkit
|
#include <itkRLEImage.h>
Run-Length Encoded image. It saves memory for label images at the expense of processing times. Unsuitable for ordinary images (in which case it is counterproductive).
It is best if pixel type and counter type have the same byte size (for memory alignment purposes).
Acknowledgement: This work is supported by NIH grant R01 EB014346, "Continued development and maintenance of the ITK-SNAP 3D image segmentation software."
Definition at line 52 of file itkRLEImage.h.
Public Types | |
typedef itk::Image< RLLine, VImageDimension-1 > | BufferType |
typedef itk::SmartPointer < const Self > | ConstPointer |
typedef itk::WeakPointer < const Self > | ConstWeakPointer |
typedef Superclass::DirectionType | DirectionType |
typedef Superclass::IndexType | IndexType |
typedef Superclass::IndexValueType | IndexValueType |
typedef RLLine | InternalPixelType |
typedef Superclass::OffsetType | OffsetType |
typedef Superclass::OffsetValueType | OffsetValueType |
typedef TPixel | PixelType |
typedef itk::SmartPointer< Self > | Pointer |
typedef Superclass::PointType | PointType |
typedef Superclass::RegionType | RegionType |
typedef CounterType | RLCounterType |
typedef std::vector< RLSegment > | RLLine |
typedef std::pair< CounterType, PixelType > | RLSegment |
typedef RLEImage | Self |
typedef Superclass::SizeType | SizeType |
typedef Superclass::SizeValueType | SizeValueType |
typedef Superclass::SpacingType | SpacingType |
typedef Superclass::SpacingValueType | SpacingValueType |
typedef itk::ImageBase < VImageDimension > | Superclass |
typedef TPixel | ValueType |
Public Types inherited from itk::ImageBase< VImageDimension > | |
typedef SmartPointer< const Self > | ConstPointer |
typedef Matrix < SpacePrecisionType, VImageDimension, VImageDimension > | DirectionType |
typedef unsigned int | ImageDimensionType |
typedef Index< VImageDimension > | IndexType |
typedef IndexType::IndexValueType | IndexValueType |
typedef Offset< VImageDimension > | OffsetType |
typedef OffsetType::OffsetValueType | OffsetValueType |
typedef SmartPointer< Self > | Pointer |
typedef Point< PointValueType, VImageDimension > | PointType |
typedef SpacePrecisionType | PointValueType |
typedef ImageRegion < VImageDimension > | RegionType |
typedef ImageBase | Self |
typedef Size< VImageDimension > | SizeType |
typedef SizeType::SizeValueType | SizeValueType |
typedef Vector < SpacingValueType, VImageDimension > | SpacingType |
typedef SpacePrecisionType | SpacingValueType |
typedef DataObject | Superclass |
Public Types inherited from itk::DataObject | |
typedef SmartPointer< const Self > | ConstPointer |
typedef std::string | DataObjectIdentifierType |
typedef std::vector< Pointer > ::size_type | DataObjectPointerArraySizeType |
typedef SmartPointer< Self > | Pointer |
typedef DataObject | Self |
typedef Object | Superclass |
Public Types inherited from itk::Object | |
typedef SmartPointer< const Self > | ConstPointer |
typedef SmartPointer< Self > | Pointer |
typedef Object | Self |
typedef LightObject | Superclass |
Public Types inherited from itk::LightObject | |
typedef SmartPointer< const Self > | ConstPointer |
typedef SmartPointer< Self > | Pointer |
typedef LightObject | Self |
Public Member Functions | |
virtual void | Allocate (bool initialize=false) |
void | CleanUp () const |
virtual ::itk::LightObject::Pointer | CreateAnother () const |
void | FillBuffer (const TPixel &value) |
BufferType::Pointer | GetBuffer () |
BufferType::Pointer | GetBuffer () const |
virtual const char * | GetNameOfClass () const |
virtual unsigned int | GetNumberOfComponentsPerPixel () const |
bool | GetOnTheFlyCleanup () const |
const TPixel & | GetPixel (const IndexType &index) const |
const TPixel & | operator[] (const IndexType &index) const |
virtual void | SetBufferedRegion (const RegionType ®ion) |
virtual void | SetLargestPossibleRegion (const RegionType ®ion) |
void | SetPixel (const IndexType &index, const TPixel &value) |
int | SetPixel (RLLine &line, IndexValueType &segmentRemainder, SizeValueType &m_RealIndex, const TPixel &value) |
virtual void | SetRequestedRegion (const RegionType ®ion) |
virtual void | Initialize () |
void | SetOnTheFlyCleanup (bool value) |
Public Member Functions inherited from itk::ImageBase< VImageDimension > | |
OffsetValueType | ComputeOffset (const IndexType &ind) const |
virtual void | CopyInformation (const DataObject *data) override |
virtual const RegionType & | GetBufferedRegion () const |
virtual const DirectionType & | GetDirection () const |
virtual const DirectionType & | GetInverseDirection () const |
virtual const RegionType & | GetLargestPossibleRegion () const |
virtual const PointType & | GetOrigin () const |
virtual const RegionType & | GetRequestedRegion () const |
virtual const SpacingType & | GetSpacing () const |
virtual void | Graft (const Self *data) |
virtual bool | RequestedRegionIsOutsideOfTheBufferedRegion () override |
virtual void | SetDirection (const DirectionType &direction) |
virtual void | SetRegions (const SizeType &size) |
virtual void | SetRequestedRegion (const DataObject *data) override |
virtual void | SetRequestedRegionToLargestPossibleRegion () override |
template<typename TCoordRep > | |
void | TransformLocalVectorToPhysicalVector (const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const |
template<typename TCoordRep , typename TIndexRep > | |
bool | TransformPhysicalPointToContinuousIndex (const Point< TCoordRep, VImageDimension > &point, ContinuousIndex< TIndexRep, VImageDimension > &index) const |
template<typename TCoordRep > | |
bool | TransformPhysicalPointToIndex (const Point< TCoordRep, VImageDimension > &point, IndexType &index) const |
template<typename TCoordRep > | |
void | TransformPhysicalVectorToLocalVector (const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const |
virtual void | UpdateOutputData () override |
virtual void | UpdateOutputInformation () override |
virtual bool | VerifyRequestedRegion () override |
virtual void | SetOrigin (PointType _arg) |
virtual void | SetOrigin (const double origin[VImageDimension]) |
virtual void | SetOrigin (const float origin[VImageDimension]) |
virtual void | SetRegions (const RegionType ®ion) |
const OffsetValueType * | GetOffsetTable () 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 , typename TIndexRep > | |
void | TransformContinuousIndexToPhysicalPoint (const ContinuousIndex< TIndexRep, VImageDimension > &index, Point< TCoordRep, VImageDimension > &point) const |
template<typename TCoordRep > | |
void | TransformIndexToPhysicalPoint (const IndexType &index, Point< TCoordRep, VImageDimension > &point) const |
virtual void | SetNumberOfComponentsPerPixel (unsigned int) |
Public Member Functions inherited from itk::DataObject | |
virtual void | DataHasBeenGenerated () |
void | DisconnectPipeline () |
bool | GetDataReleased () const |
virtual const bool & | GetReleaseDataFlag () const |
SmartPointerForwardReference < ProcessObject > | GetSource () const |
DataObjectPointerArraySizeType | GetSourceOutputIndex () const |
const DataObjectIdentifierType & | GetSourceOutputName () const |
virtual ModifiedTimeType | GetUpdateMTime () const |
virtual void | PrepareForNewData () |
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 ModifiedTimeType & | GetPipelineMTime () const |
virtual void | SetRealTimeStamp (RealTimeStamp _arg) |
virtual const RealTimeStamp & | GetRealTimeStamp () const |
Public Member Functions inherited from itk::Object | |
unsigned long | AddObserver (const EventObject &event, Command *) |
unsigned long | AddObserver (const EventObject &event, Command *) const |
virtual void | DebugOff () const |
virtual void | DebugOn () const |
Command * | GetCommand (unsigned long tag) |
bool | GetDebug () const |
MetaDataDictionary & | GetMetaDataDictionary () |
const MetaDataDictionary & | GetMetaDataDictionary () const |
virtual ModifiedTimeType | GetMTime () const |
virtual const TimeStamp & | GetTimeStamp () const |
bool | HasObserver (const EventObject &event) const |
void | InvokeEvent (const EventObject &) |
void | InvokeEvent (const EventObject &) const |
virtual void | Modified () const |
virtual void | Register () const override |
void | RemoveAllObservers () |
void | RemoveObserver (unsigned long tag) |
void | SetDebug (bool debugFlag) const |
void | SetMetaDataDictionary (const MetaDataDictionary &rhs) |
virtual void | SetReferenceCount (int) override |
virtual void | UnRegister () const noexceptoverride |
virtual void | SetObjectName (std::string _arg) |
virtual const std::string & | GetObjectName () const |
Public Member Functions inherited from itk::LightObject | |
virtual void | Delete () |
virtual int | GetReferenceCount () const |
itkCloneMacro (Self) | |
void | Print (std::ostream &os, Indent indent=0) const |
Static Public Member Functions | |
static Pointer | New () |
static BufferType::IndexType | truncateIndex (const IndexType &index) |
static BufferType::RegionType | truncateRegion (const RegionType ®ion) |
static BufferType::SizeType | truncateSize (const SizeType &size) |
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 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 flag) |
Static Public Member Functions inherited from itk::LightObject | |
static void | BreakOnError () |
static Pointer | New () |
Static Public Attributes | |
static const unsigned int | ImageDimension = VImageDimension |
Static Public Attributes inherited from itk::ImageBase< VImageDimension > | |
static const ImageDimensionType | ImageDimension = VImageDimension |
Protected Member Functions | |
void | CleanUpLine (RLLine &line) const |
virtual void | ComputeIndexToPhysicalPointMatrices () |
void | PrintSelf (std::ostream &os, itk::Indent indent) const |
RLEImage () | |
virtual | ~RLEImage () |
Protected Member Functions inherited from itk::ImageBase< VImageDimension > | |
void | ComputeOffsetTable () |
virtual void | Graft (const DataObject *data) override |
ImageBase () | |
virtual void | InitializeBufferedRegion () |
~ImageBase () override | |
OffsetValueType | FastComputeOffset (const IndexType &ind) const |
IndexType | FastComputeIndex (OffsetValueType offset) const |
Protected Member Functions inherited from itk::DataObject | |
DataObject () | |
virtual void | PropagateResetPipeline () |
virtual | ~DataObject () override |
Protected Member Functions inherited from itk::Object | |
Object () | |
bool | PrintObservers (std::ostream &os, Indent indent) const |
virtual void | SetTimeStamp (const TimeStamp &time) |
virtual | ~Object () override |
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 () |
Private Member Functions | |
void | operator= (const Self &) |
RLEImage (const Self &) | |
Private Attributes | |
BufferType::Pointer | m_Buffer |
bool | m_OnTheFlyCleanup |
Additional Inherited Members | |
Protected Attributes inherited from itk::ImageBase< VImageDimension > | |
DirectionType | m_Direction |
DirectionType | m_IndexToPhysicalPoint |
DirectionType | m_InverseDirection |
PointType | m_Origin |
DirectionType | m_PhysicalPointToIndex |
SpacingType | m_Spacing |
Protected Attributes inherited from itk::LightObject | |
AtomicInt< int > | m_ReferenceCount |
typedef itk::Image< RLLine, VImageDimension - 1 > itk::RLEImage< TPixel, VImageDimension, CounterType >::BufferType |
Typedef for the internally used buffer.
Definition at line 222 of file itkRLEImage.h.
typedef itk::SmartPointer< const Self > itk::RLEImage< TPixel, VImageDimension, CounterType >::ConstPointer |
Definition at line 60 of file itkRLEImage.h.
typedef itk::WeakPointer< const Self > itk::RLEImage< TPixel, VImageDimension, CounterType >::ConstWeakPointer |
Definition at line 61 of file itkRLEImage.h.
typedef Superclass::DirectionType itk::RLEImage< TPixel, VImageDimension, CounterType >::DirectionType |
Direction typedef support. A matrix of direction cosines.
Definition at line 111 of file itkRLEImage.h.
typedef Superclass::IndexType itk::RLEImage< TPixel, VImageDimension, CounterType >::IndexType |
Index typedef support. An index is used to access pixel values.
Definition at line 100 of file itkRLEImage.h.
typedef Superclass::IndexValueType itk::RLEImage< TPixel, VImageDimension, CounterType >::IndexValueType |
Definition at line 101 of file itkRLEImage.h.
typedef RLLine itk::RLEImage< TPixel, VImageDimension, CounterType >::InternalPixelType |
Internal Pixel representation. Used to maintain a uniform API with Image Adaptors and allow to keep a particular internal representation of data while showing a different external representation.
Definition at line 89 of file itkRLEImage.h.
typedef Superclass::OffsetType itk::RLEImage< TPixel, VImageDimension, CounterType >::OffsetType |
Offset typedef support. An offset is used to access pixel values.
Definition at line 104 of file itkRLEImage.h.
typedef Superclass::OffsetValueType itk::RLEImage< TPixel, VImageDimension, CounterType >::OffsetValueType |
Offset typedef (relative position between indices)
Definition at line 127 of file itkRLEImage.h.
typedef TPixel itk::RLEImage< TPixel, VImageDimension, CounterType >::PixelType |
Pixel typedef support. Used to declare pixel type in filters or other operations.
Definition at line 67 of file itkRLEImage.h.
typedef itk::SmartPointer< Self > itk::RLEImage< TPixel, VImageDimension, CounterType >::Pointer |
Definition at line 59 of file itkRLEImage.h.
typedef Superclass::PointType itk::RLEImage< TPixel, VImageDimension, CounterType >::PointType |
Origin typedef support. The origin is the geometric coordinates of the index (0,0).
Definition at line 124 of file itkRLEImage.h.
typedef Superclass::RegionType itk::RLEImage< TPixel, VImageDimension, CounterType >::RegionType |
Region typedef support. A region is used to specify a subset of an image.
Definition at line 115 of file itkRLEImage.h.
typedef CounterType itk::RLEImage< TPixel, VImageDimension, CounterType >::RLCounterType |
Definition at line 73 of file itkRLEImage.h.
typedef std::vector< RLSegment > itk::RLEImage< TPixel, VImageDimension, CounterType >::RLLine |
A Run-Length encoded line of pixels.
Definition at line 83 of file itkRLEImage.h.
typedef std::pair< CounterType, PixelType > itk::RLEImage< TPixel, VImageDimension, CounterType >::RLSegment |
First element is count of repetitions, second element is the pixel value.
Definition at line 80 of file itkRLEImage.h.
typedef RLEImage itk::RLEImage< TPixel, VImageDimension, CounterType >::Self |
Standard class typedefs
Definition at line 57 of file itkRLEImage.h.
typedef Superclass::SizeType itk::RLEImage< TPixel, VImageDimension, CounterType >::SizeType |
Size typedef support. A size is used to define region bounds.
Definition at line 107 of file itkRLEImage.h.
typedef Superclass::SizeValueType itk::RLEImage< TPixel, VImageDimension, CounterType >::SizeValueType |
Definition at line 108 of file itkRLEImage.h.
typedef Superclass::SpacingType itk::RLEImage< TPixel, VImageDimension, CounterType >::SpacingType |
Spacing typedef support. Spacing holds the size of a pixel. The spacing is the geometric distance between image samples.
Definition at line 119 of file itkRLEImage.h.
typedef Superclass::SpacingValueType itk::RLEImage< TPixel, VImageDimension, CounterType >::SpacingValueType |
Definition at line 120 of file itkRLEImage.h.
typedef itk::ImageBase< VImageDimension > itk::RLEImage< TPixel, VImageDimension, CounterType >::Superclass |
Definition at line 58 of file itkRLEImage.h.
typedef TPixel itk::RLEImage< TPixel, VImageDimension, CounterType >::ValueType |
Typedef alias for PixelType
Definition at line 76 of file itkRLEImage.h.
|
inlineprotected |
Definition at line 281 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::m_Buffer, and itk::Image< TPixel, VImageDimension >::New().
|
inlineprotectedvirtual |
Definition at line 292 of file itkRLEImage.h.
|
private |
|
virtual |
Allocate the image memory. The size of the image must already be set, e.g. by calling SetRegions(). Pixel values are initialized using default constructor.
Reimplemented from itk::ImageBase< VImageDimension >.
void itk::RLEImage< TPixel, VImageDimension, CounterType >::CleanUp | ( | ) | const |
Merges adjacent segments with duplicate values. Automatically called when turning on OnTheFlyCleanup.
Referenced by itk::RLEImage< TPixel, VImageDimension, CounterType >::SetOnTheFlyCleanup().
|
protected |
Merges adjacent segments with duplicate values in a single line.
|
inlineprotectedvirtual |
Compute helper matrices used to transform Index coordinates to PhysicalPoint coordinates and back. This method is virtual and will be overloaded in derived classes in order to provide backward compatibility behavior in classes that did not used to take image orientation into account.
Reimplemented from itk::ImageBase< VImageDimension >.
Definition at line 299 of file itkRLEImage.h.
References itk::ImageBase< VImageDimension >::ComputeIndexToPhysicalPointMatrices().
|
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 >.
void itk::RLEImage< TPixel, VImageDimension, CounterType >::FillBuffer | ( | const TPixel & | value | ) |
Fill the image buffer with a value. Be sure to call Allocate() first.
|
inline |
We need to allow itk-style iterators to be constructed.
Definition at line 226 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::m_Buffer.
|
inline |
We need to allow itk-style const iterators to be constructed.
Definition at line 233 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::m_Buffer.
|
virtual |
Run-time type information (and related methods).
Reimplemented from itk::ImageBase< VImageDimension >.
|
inlinevirtual |
INTERNAL This method is used internally by filters to copy meta-data from the output to the input. Users should not have a need to use this method.
Filters that override the ProcessObject's GenerateOutputInformation() should generally have the following line if they want to propagate meta- data for both Image and VectorImage
Returns/Sets the number of components in the image. Note that in the ImageBase implementation, this always returns 1. Image returns the
returns the vector length set by the user.
Reimplemented from itk::ImageBase< VImageDimension >.
Definition at line 212 of file itkRLEImage.h.
References itk::NumericTraits< T >::GetLength().
|
inline |
Should same-valued segments be merged on the fly? On the fly merging usually provides better performance.
Definition at line 258 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::m_OnTheFlyCleanup.
const TPixel& itk::RLEImage< TPixel, VImageDimension, CounterType >::GetPixel | ( | const IndexType & | index | ) | const |
Get a pixel. SLOW! Better use iterators for pixel access.
Referenced by itk::RLEImage< TPixel, VImageDimension, CounterType >::operator[]().
|
inlinevirtual |
Restore the data object to its initial state. This means releasing memory.
Reimplemented from itk::ImageBase< VImageDimension >.
Definition at line 138 of file itkRLEImage.h.
References itk::ImageBase< VImageDimension >::Initialize(), itk::RLEImage< TPixel, VImageDimension, CounterType >::m_Buffer, itk::RLEImage< TPixel, VImageDimension, CounterType >::m_OnTheFlyCleanup, and itk::Image< TPixel, VImageDimension >::New().
|
static |
Method for creation through the object factory.
|
private |
|
inline |
Access a pixel. Chaning it changes the whole RLE segment!
Get a reference to a pixel. Chaning it changes the whole RLE segment!Access a pixel. This version can only be an rvalue. SLOW -> Use iterators instead.
Definition at line 206 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::GetPixel().
|
protectedvirtual |
Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.
Reimplemented from itk::ImageBase< VImageDimension >.
|
inlinevirtual |
Set the region object that defines the size and starting index of the region of the image currently loaded in memory.
Reimplemented from itk::ImageBase< VImageDimension >.
Definition at line 160 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::m_Buffer, itk::ImageBase< VImageDimension >::SetBufferedRegion(), and itk::RLEImage< TPixel, VImageDimension, CounterType >::truncateRegion().
|
inlinevirtual |
Set the region object that defines the size and starting index for the largest possible region this image could represent. This is used in determining how much memory would be needed to load an entire dataset. It is also used to determine boundary true conditions.
Reimplemented from itk::ImageBase< VImageDimension >.
Definition at line 153 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::m_Buffer, itk::ImageBase< VImageDimension >::SetLargestPossibleRegion(), and itk::RLEImage< TPixel, VImageDimension, CounterType >::truncateRegion().
|
inline |
Should same-valued segments be merged on the fly? On the fly merging usually provides better performance.
Definition at line 266 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::CleanUp(), and itk::RLEImage< TPixel, VImageDimension, CounterType >::m_OnTheFlyCleanup.
void itk::RLEImage< TPixel, VImageDimension, CounterType >::SetPixel | ( | const IndexType & | index, |
const TPixel & | value | ||
) |
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. SLOW -> Use iterators instead.
int itk::RLEImage< TPixel, VImageDimension, CounterType >::SetPixel | ( | RLLine & | line, |
IndexValueType & | segmentRemainder, | ||
SizeValueType & | m_RealIndex, | ||
const TPixel & | value | ||
) |
Set a pixel value in the given line and updates segmentRemainder and m_RealIndex to refer to the same pixel. Returns difference in line length which happens due to merging or splitting segments. This method is used by iterators directly.
|
inlinevirtual |
Set the region object that defines the size and starting index for the region of the image requested (i.e., the region of the image to be operated on by a filter). Setting the RequestedRegion does not cause the object to be modified. This method is called internally by the pipeline and therefore bypasses the modified time calculation.
Reimplemented from itk::ImageBase< VImageDimension >.
Definition at line 169 of file itkRLEImage.h.
References itk::RLEImage< TPixel, VImageDimension, CounterType >::m_Buffer, itk::ImageBase< VImageDimension >::SetRequestedRegion(), and itk::RLEImage< TPixel, VImageDimension, CounterType >::truncateRegion().
|
inlinestatic |
Returns N-1-dimensional index, the remainder after 0-index is removed.
|
static |
Returns N-1-dimensional region, the remainder after 0-index and size are removed.
Referenced by itk::RLEImage< TPixel, VImageDimension, CounterType >::SetBufferedRegion(), itk::RLEImage< TPixel, VImageDimension, CounterType >::SetLargestPossibleRegion(), and itk::RLEImage< TPixel, VImageDimension, CounterType >::SetRequestedRegion().
|
inlinestatic |
Returns N-1-dimensional size, the remainder after 0-size is removed.
|
static |
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 97 of file itkRLEImage.h.
|
mutableprivate |
Memory for the current buffer.
Definition at line 316 of file itkRLEImage.h.
Referenced by itk::RLEImage< TPixel, VImageDimension, CounterType >::GetBuffer(), itk::RLEImage< TPixel, VImageDimension, CounterType >::Initialize(), itk::RLEImage< TPixel, VImageDimension, CounterType >::RLEImage(), itk::RLEImage< TPixel, VImageDimension, CounterType >::SetBufferedRegion(), itk::RLEImage< TPixel, VImageDimension, CounterType >::SetLargestPossibleRegion(), and itk::RLEImage< TPixel, VImageDimension, CounterType >::SetRequestedRegion().
|
private |