ITK  6.0.0
Insight Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | Friends | List of all members
itk::ImageRegion< VImageDimension > Class Template Referencefinal

#include <itkImageRegion.h>

Detailed Description

template<unsigned int VImageDimension>
class itk::ImageRegion< VImageDimension >

An image region represents a structured region of data.

ImageRegion is an class that represents some structured portion or piece of an Image. The ImageRegion is represented with an index and a size in each of the n-dimensions of the image. (The index is the corner of the image, the size is the lengths of the image in each of the topological directions.)

See also
Region
Index
Size
MeshRegion
ITK Sphinx Examples:
Examples
SphinxExamples/src/Core/Common/CreateAnImageRegion/Code.cxx, SphinxExamples/src/Core/Common/ImageBufferAndIndexRange/Code.cxx, SphinxExamples/src/Core/Common/IsPixelInsideRegion/Code.cxx, SphinxExamples/src/Core/Common/IterateImageStartingAtSeed/Code.cxx, SphinxExamples/src/Core/Common/NeighborhoodIteratorOnVectorImage/Code.cxx, SphinxExamples/src/Core/Common/PassImageToFunction/Code.cxx, SphinxExamples/src/Core/Common/ReturnObjectFromFunction/Code.cxx, SphinxExamples/src/Core/ImageFunction/ComputeMedianOfImageAtPixel/Code.cxx, SphinxExamples/src/Core/Transform/TranslateAVectorImage/Code.cxx, SphinxExamples/src/Filtering/BinaryMathematicalMorphology/ClosingBinaryImage/Code.cxx, SphinxExamples/src/Filtering/BinaryMathematicalMorphology/OpeningBinaryImage/Code.cxx, SphinxExamples/src/Filtering/DistanceMap/ApproxDistanceMapOfBinary/Code.cxx, SphinxExamples/src/Filtering/DistanceMap/MaurerDistanceMapOfBinary/Code.cxx, SphinxExamples/src/Filtering/DistanceMap/MeanDistanceBetweenAllPointsOnTwoCurves/Code.cxx, SphinxExamples/src/Filtering/DistanceMap/SignedDistanceMapOfBinary/Code.cxx, SphinxExamples/src/Filtering/FFT/ComputeInverseFFTOfImage/Code.cxx, SphinxExamples/src/Filtering/ImageFeature/ApplyAFilterOnlyToASpecifiedImageRegion/Code.cxx, SphinxExamples/src/Filtering/ImageFeature/FindZeroCrossingsInSignedImage/Code.cxx, SphinxExamples/src/Filtering/ImageFilterBase/ApplyKernelToEveryPixelInNonZeroImage/Code.cxx, SphinxExamples/src/Filtering/ImageGrid/RunImageFilterOnRegionOfImage/Code.cxx, SphinxExamples/src/Filtering/Path/ExtractContoursFromImage/Code.cxx, SphinxExamples/src/ImageCompareCommand.cxx, SphinxExamples/src/IO/ImageBase/ConvertImageToAnotherType/Code.cxx, SphinxExamples/src/Numerics/Statistics/ComputeHistogramOfMaskedRegion/Code.cxx, SphinxExamples/src/Numerics/Statistics/ComputeTextureFeatures/Code.cxx, and SphinxExamples/src/Segmentation/Watersheds/MorphologicalWatershedSegmentation/Code.cxx.

Definition at line 80 of file itkImageRegion.h.

+ Collaboration diagram for itk::ImageRegion< VImageDimension >:

Public Types

using IndexType = Index< VImageDimension >
 
using IndexValueArrayType = IndexValueType[ImageDimension]
 
using IndexValueType = typename IndexType::IndexValueType
 
using OffsetTableType = OffsetValueType[ImageDimension+1]
 
using OffsetType = typename IndexType::OffsetType
 
using OffsetValueType = typename OffsetType::OffsetValueType
 
using Self = ImageRegion
 
using SizeType = Size< VImageDimension >
 
using SizeValueType = typename SizeType::SizeValueType
 
using SliceRegion = ImageRegion< Self::SliceDimension >
 

Public Member Functions

void ComputeOffsetTable (OffsetTableType offsetTable) const
 
bool Crop (const Self &region)
 
const char * GetNameOfClass () const itkRegionOverrideMacro
 
SizeValueType GetNumberOfPixels () const
 
Region::RegionEnum GetRegionType () const itkRegionOverrideMacro
 
IndexType GetUpperIndex () const
 
 ImageRegion () noexcept=default
 
 ImageRegion (const IndexType &index, const SizeType &size) noexcept
 
 ImageRegion (const Self &) noexcept=default
 
 ImageRegion (const SizeType &size) noexcept
 
 ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION (Self)
 
Selfoperator= (const Self &) noexcept=default
 
bool operator== (const Self &region) const noexcept
 
void PadByRadius (const IndexValueArrayType radius)
 
void PadByRadius (const SizeType &radius)
 
void PadByRadius (OffsetValueType radius)
 
void Print (std::ostream &os, Indent indent=0) const itkRegionOverrideMacro
 
void SetIndex (const IndexType &index)
 
void SetSize (const SizeType &size)
 
void SetUpperIndex (const IndexType &idx)
 
bool ShrinkByRadius (const IndexValueArrayType radius)
 
bool ShrinkByRadius (const SizeType &radius)
 
bool ShrinkByRadius (OffsetValueType radius)
 
SliceRegion Slice (const unsigned int dim) const
 
 ~ImageRegion () itkRegionOverrideMacro=default
 
const IndexTypeGetIndex () const
 
IndexTypeGetModifiableIndex ()
 
const SizeTypeGetSize () const
 
SizeTypeGetModifiableSize ()
 
void SetSize (unsigned int i, SizeValueType sze)
 
SizeValueType GetSize (unsigned int i) const
 
void SetIndex (unsigned int i, IndexValueType sze)
 
IndexValueType GetIndex (unsigned int i) const
 
bool IsInside (const IndexType &index) const
 
template<typename TCoordinateType >
bool IsInside (const ContinuousIndex< TCoordinateType, VImageDimension > &index) const
 
bool IsInside (const Self &otherRegion) const
 
template<vcl_size_t VTupleIndex>
auto & get ()
 
template<vcl_size_t VTupleIndex>
const auto & get () const
 

Static Public Member Functions

static constexpr unsigned int GetImageDimension ()
 

Static Public Attributes

static constexpr unsigned int ImageDimension = VImageDimension
 
static constexpr unsigned int SliceDimension = ImageDimension - (ImageDimension > 1)
 

Protected Member Functions

void PrintSelf (std::ostream &os, Indent indent) const itkRegionOverrideMacro
 

Private Attributes

IndexType m_Index = { { 0 } }
 
SizeType m_Size = { { 0 } }
 

Friends

class ImageBase< VImageDimension >
 

Member Typedef Documentation

◆ IndexType

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::IndexType = Index<VImageDimension>

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

Definition at line 116 of file itkImageRegion.h.

◆ IndexValueArrayType

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::IndexValueArrayType = IndexValueType[ImageDimension]

Definition at line 120 of file itkImageRegion.h.

◆ IndexValueType

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::IndexValueType = typename IndexType::IndexValueType

Definition at line 117 of file itkImageRegion.h.

◆ OffsetTableType

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::OffsetTableType = OffsetValueType[ImageDimension + 1]

Definition at line 121 of file itkImageRegion.h.

◆ OffsetType

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::OffsetType = typename IndexType::OffsetType

Definition at line 118 of file itkImageRegion.h.

◆ OffsetValueType

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::OffsetValueType = typename OffsetType::OffsetValueType

Definition at line 119 of file itkImageRegion.h.

◆ Self

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::Self = ImageRegion

Standard class type aliases.

Definition at line 88 of file itkImageRegion.h.

◆ SizeType

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::SizeType = Size<VImageDimension>

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

Definition at line 124 of file itkImageRegion.h.

◆ SizeValueType

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::SizeValueType = typename SizeType::SizeValueType

Definition at line 125 of file itkImageRegion.h.

◆ SliceRegion

template<unsigned int VImageDimension>
using itk::ImageRegion< VImageDimension >::SliceRegion = ImageRegion<Self::SliceDimension>

Slice region type alias. SliceRegion is one dimension less than Self.

Definition at line 128 of file itkImageRegion.h.

Constructor & Destructor Documentation

◆ ImageRegion() [1/4]

template<unsigned int VImageDimension>
itk::ImageRegion< VImageDimension >::ImageRegion ( )
defaultnoexcept

Constructor. ImageRegion is a lightweight object that is not reference counted, so the constructor is public. Its two data members are filled with zeros (using C++11 default member initializers).

◆ ~ImageRegion()

template<unsigned int VImageDimension>
itk::ImageRegion< VImageDimension >::~ImageRegion ( )
default

Destructor. ImageRegion is a lightweight object that is not reference counted, so the destructor is public.

◆ ImageRegion() [2/4]

template<unsigned int VImageDimension>
itk::ImageRegion< VImageDimension >::ImageRegion ( const Self )
defaultnoexcept

Copy constructor. ImageRegion is a lightweight object that is not reference counted, so the copy constructor is public.

◆ ImageRegion() [3/4]

template<unsigned int VImageDimension>
itk::ImageRegion< VImageDimension >::ImageRegion ( const IndexType index,
const SizeType size 
)
inlinenoexcept

Constructor that takes an index and size. ImageRegion is a lightweight object that is not reference counted, so this constructor is public.

Note
This constructor supports class template argument deduction (CTAD).

Definition at line 157 of file itkImageRegion.h.

◆ ImageRegion() [4/4]

template<unsigned int VImageDimension>
itk::ImageRegion< VImageDimension >::ImageRegion ( const SizeType size)
inlinenoexcept

Constructor that takes a size and assumes an index of zeros. ImageRegion is lightweight object that is not reference counted so this constructor is public.

Note
This constructor supports class template argument deduction (CTAD).

Definition at line 168 of file itkImageRegion.h.

Member Function Documentation

◆ ComputeOffsetTable()

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::ComputeOffsetTable ( OffsetTableType  offsetTable) const

Compute an offset table based on the Size.

◆ Crop()

template<unsigned int VImageDimension>
bool itk::ImageRegion< VImageDimension >::Crop ( const Self region)

Crop a region by another region. If this region is outside of the crop, this method returns false and does not modify the region. Otherwise, this method returns true and the region is modified to reflect the crop.

Examples
SphinxExamples/src/Core/Common/ImageRegionOverlap/Code.cxx.

◆ get() [1/2]

template<unsigned int VImageDimension>
template<vcl_size_t VTupleIndex>
auto& itk::ImageRegion< VImageDimension >::get ( )
inline

Supports tuple-like access: get<0>() returns a reference to the index and get<1>() returns a reference to the size of the region.

Definition at line 375 of file itkImageRegion.h.

◆ get() [2/2]

template<unsigned int VImageDimension>
template<vcl_size_t VTupleIndex>
const auto& itk::ImageRegion< VImageDimension >::get ( ) const
inline

Supports tuple-like access. Const overload.

Definition at line 392 of file itkImageRegion.h.

◆ GetImageDimension()

template<unsigned int VImageDimension>
static constexpr unsigned int itk::ImageRegion< VImageDimension >::GetImageDimension ( )
inlinestaticconstexpr

Dimension of the image available at compile-time and at run time.

Definition at line 110 of file itkImageRegion.h.

◆ GetIndex() [1/2]

template<unsigned int VImageDimension>
const IndexType& itk::ImageRegion< VImageDimension >::GetIndex ( ) const
inline

◆ GetIndex() [2/2]

template<unsigned int VImageDimension>
IndexValueType itk::ImageRegion< VImageDimension >::GetIndex ( unsigned int  i) const
inline

Convenience methods to get and set the index of the particular dimension i.

Definition at line 242 of file itkImageRegion.h.

◆ GetModifiableIndex()

template<unsigned int VImageDimension>
IndexType& itk::ImageRegion< VImageDimension >::GetModifiableIndex ( )
inline

Get index defining the corner of the region.

Definition at line 193 of file itkImageRegion.h.

Referenced by itk::ImageRegionSplitterBase::GetSplit().

◆ GetModifiableSize()

template<unsigned int VImageDimension>
SizeType& itk::ImageRegion< VImageDimension >::GetModifiableSize ( )
inline

Get the size of the region.

Definition at line 214 of file itkImageRegion.h.

Referenced by itk::ImageRegionSplitterBase::GetSplit().

◆ GetNameOfClass()

template<unsigned int VImageDimension>
const char* itk::ImageRegion< VImageDimension >::GetNameOfClass ( ) const
inline

Standard part of all itk objects.

Definition at line 96 of file itkImageRegion.h.

◆ GetNumberOfPixels()

template<unsigned int VImageDimension>
SizeValueType itk::ImageRegion< VImageDimension >::GetNumberOfPixels ( ) const

◆ GetRegionType()

template<unsigned int VImageDimension>
Region::RegionEnum itk::ImageRegion< VImageDimension >::GetRegionType ( ) const
inline

Return the region type. Images are described with structured regions.

Definition at line 132 of file itkImageRegion.h.

◆ GetSize() [1/2]

template<unsigned int VImageDimension>
const SizeType& itk::ImageRegion< VImageDimension >::GetSize ( ) const
inline

◆ GetSize() [2/2]

template<unsigned int VImageDimension>
SizeValueType itk::ImageRegion< VImageDimension >::GetSize ( unsigned int  i) const
inline

Convenience methods to get and set the size of the particular dimension i.

Definition at line 228 of file itkImageRegion.h.

◆ GetUpperIndex()

template<unsigned int VImageDimension>
IndexType itk::ImageRegion< VImageDimension >::GetUpperIndex ( ) const

Get index defining the upper corner of the region.

◆ IsInside() [1/3]

template<unsigned int VImageDimension>
template<typename TCoordinateType >
bool itk::ImageRegion< VImageDimension >::IsInside ( const ContinuousIndex< TCoordinateType, VImageDimension > &  index) const
inline

Test if a continuous index is inside the region. We take into account the fact that each voxel has its center at the integer coordinate and extends half way to the next integer coordinate, inclusive on all sides.

Definition at line 290 of file itkImageRegion.h.

◆ IsInside() [2/3]

template<unsigned int VImageDimension>
bool itk::ImageRegion< VImageDimension >::IsInside ( const IndexType index) const
inline

Test if an index is inside

Definition at line 271 of file itkImageRegion.h.

◆ IsInside() [3/3]

template<unsigned int VImageDimension>
bool itk::ImageRegion< VImageDimension >::IsInside ( const Self otherRegion) const
inline

Test if a region (the argument) is completely inside of this region. If the region that is passed as argument to this method, has a size of value zero, then it will not be considered to be inside of the current region, even its starting index is inside.

Definition at line 310 of file itkImageRegion.h.

◆ ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION()

template<unsigned int VImageDimension>
itk::ImageRegion< VImageDimension >::ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION ( Self  )

◆ operator=()

template<unsigned int VImageDimension>
Self& itk::ImageRegion< VImageDimension >::operator= ( const Self )
defaultnoexcept

operator=. ImageRegion is a lightweight object that is not reference counted, so operator= is public.

◆ operator==()

template<unsigned int VImageDimension>
bool itk::ImageRegion< VImageDimension >::operator== ( const Self region) const
inlinenoexcept

Compare two regions.

Definition at line 262 of file itkImageRegion.h.

◆ PadByRadius() [1/3]

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::PadByRadius ( const IndexValueArrayType  radius)

◆ PadByRadius() [2/3]

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::PadByRadius ( const SizeType radius)

◆ PadByRadius() [3/3]

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::PadByRadius ( OffsetValueType  radius)

Pad an image region by the specified radius. Region can be padded uniformly in all dimensions or can be padded by different amounts in each dimension.

◆ Print()

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::Print ( std::ostream &  os,
Indent  indent = 0 
) const

◆ PrintSelf()

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
protected

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.

◆ SetIndex() [1/2]

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::SetIndex ( const IndexType index)
inline

◆ SetIndex() [2/2]

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::SetIndex ( unsigned int  i,
IndexValueType  sze 
)
inline

Convenience methods to get and set the index of the particular dimension i.

Definition at line 237 of file itkImageRegion.h.

◆ SetSize() [1/2]

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::SetSize ( const SizeType size)
inline

Set the size of the region. This plus the index determines the rectangular shape, or extent, of the region.

Examples
Examples/DataRepresentation/Image/Image1.cxx, Examples/DataRepresentation/Image/Image3.cxx, Examples/DataRepresentation/Image/Image4.cxx, Examples/DataRepresentation/Image/VectorImage.cxx, Examples/Filtering/DigitallyReconstructedRadiograph1.cxx, Examples/IO/ImageReadExtractFilterInsertWrite.cxx, Examples/IO/ImageReadExtractWrite.cxx, Examples/IO/ImageReadRegionOfInterestWrite.cxx, Examples/RegistrationITKv4/ImageRegistration8.cxx, Examples/Segmentation/GibbsPriorImageFilter1.cxx, Examples/Segmentation/HoughTransform2DCirclesImageFilter.cxx, Examples/SpatialObjects/ImageMaskSpatialObject.cxx, Examples/SpatialObjects/ImageSpatialObject.cxx, SphinxExamples/src/Core/Common/CreateAnImage/Code.cxx, SphinxExamples/src/Core/Common/CreateAnImageOfVectors/Code.cxx, SphinxExamples/src/Core/Common/IterateLineThroughImage/Code.cxx, SphinxExamples/src/Core/Common/IterateLineThroughImageWithoutWriteAccess/Code.cxx, SphinxExamples/src/Core/Common/IterateRegionWithAccessToIndexWithoutWriteAccess/Code.cxx, SphinxExamples/src/Core/Common/IterateRegionWithAccessToIndexWithWriteAccess/Code.cxx, SphinxExamples/src/Core/Common/IterateRegionWithNeighborhood/Code.cxx, SphinxExamples/src/Core/Common/IterateRegionWithoutWriteAccess/Code.cxx, SphinxExamples/src/Core/Common/IterateRegionWithWriteAccess/Code.cxx, SphinxExamples/src/Core/Common/IterateWithNeighborhoodWithoutAccess/Code.cxx, SphinxExamples/src/Core/Common/MakeOutOfBoundsPixelsReturnConstValue/Code.cxx, SphinxExamples/src/Core/Common/RandomSelectOfPixelsFromRegion/Code.cxx, SphinxExamples/src/Core/Common/StoreNonPixelDataInImage/Code.cxx, SphinxExamples/src/Core/ImageAdaptors/AddConstantToPixelsWithoutDuplicatingImage/Code.cxx, SphinxExamples/src/Core/ImageAdaptors/ExtractChannelOfImageWithMultipleComponents/Code.cxx, SphinxExamples/src/Core/ImageAdaptors/PresentImageAfterOperation/Code.cxx, SphinxExamples/src/Core/ImageAdaptors/ProcessNthComponentOfVectorImage/Code.cxx, SphinxExamples/src/Core/ImageAdaptors/ViewComponentVectorImageAsScaleImage/Code.cxx, SphinxExamples/src/Filtering/Convolution/ConvolveImageWithKernel/Code.cxx, SphinxExamples/src/Filtering/FastMarching/CreateDistanceMapFromSeeds/Code.cxx, SphinxExamples/src/Filtering/ImageCompare/AbsValueOfTwoImages/Code.cxx, SphinxExamples/src/Filtering/ImageCompare/CombineTwoImagesWithCheckerBoardPattern/Code.cxx, SphinxExamples/src/Filtering/ImageCompare/SquaredDifferenceOfTwoImages/Code.cxx, SphinxExamples/src/Filtering/ImageCompose/ComposeVectorFromThreeScalarImages/Code.cxx, SphinxExamples/src/Filtering/ImageFusion/ColorBoundariesOfRegions/Code.cxx, SphinxExamples/src/Filtering/ImageFusion/ColorLabeledRegions/Code.cxx, SphinxExamples/src/Filtering/ImageFusion/OverlayLabelMapOnImage/Code.cxx, SphinxExamples/src/Filtering/ImageGradient/GradientOfVectorImage/Code.cxx, SphinxExamples/src/Filtering/ImageGrid/ProcessA2DSliceOfA3DImage/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/AddConstantToEveryPixel/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/BinaryANDTwoImages/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/BinaryORTwoImages/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/BinaryXORTwoImages/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/CompareTwoImagesAndSetOutputPixelToMax/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/CompareTwoImagesAndSetOutputPixelToMin/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/SquareEveryPixel/Code.cxx, SphinxExamples/src/Filtering/ImageIntensity/SubtractConstantFromEveryPixel/Code.cxx, SphinxExamples/src/Filtering/ImageStatistics/StatisticalPropertiesOfRegions/Code.cxx, SphinxExamples/src/Filtering/LabelMap/ConvertLabelMapToImage/Code.cxx, SphinxExamples/src/Filtering/LabelMap/InvertImageUsingBinaryNot/Code.cxx, SphinxExamples/src/Filtering/LabelMap/KeepRegionsAboveLevel/Code.cxx, SphinxExamples/src/Filtering/LabelMap/KeepRegionsThatMeetSpecific/Code.cxx, SphinxExamples/src/Filtering/LabelMap/LabelBinaryRegionsAndGetProperties/Code.cxx, SphinxExamples/src/Filtering/LabelMap/LabelBinaryRegionsInImage/Code.cxx, SphinxExamples/src/Filtering/Thresholding/SeparateGroundUsingOtsu/Code.cxx, SphinxExamples/src/Nonunit/Review/GeometricPropertiesOfRegion/Code.cxx, and SphinxExamples/src/Registration/Common/RegisterImageToAnotherUsingLandmarks/Code.cxx.

Definition at line 202 of file itkImageRegion.h.

Referenced by itk::BoxMeanCalculatorFunction(), itk::BoxSigmaCalculatorFunction(), itk::ImageIORegionAdaptor< VDimension >::Convert(), itk::ImageToImageFilterDetail::ExtractImageFilterCopyRegion(), itk::ImageToImageFilterDetail::ImageToImageFilterDefaultCopyRegion(), itk::MultiThreaderBase::ParallelizeImageRegion(), and itk::ImageBase< TImage::ImageDimension >::SetRegions().

◆ SetSize() [2/2]

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::SetSize ( unsigned int  i,
SizeValueType  sze 
)
inline

Convenience methods to get and set the size of the particular dimension i.

Definition at line 223 of file itkImageRegion.h.

◆ SetUpperIndex()

template<unsigned int VImageDimension>
void itk::ImageRegion< VImageDimension >::SetUpperIndex ( const IndexType idx)

Modify the Size of the ImageRegion so that the provided index will be the upper corner index.

◆ ShrinkByRadius() [1/3]

template<unsigned int VImageDimension>
bool itk::ImageRegion< VImageDimension >::ShrinkByRadius ( const IndexValueArrayType  radius)

◆ ShrinkByRadius() [2/3]

template<unsigned int VImageDimension>
bool itk::ImageRegion< VImageDimension >::ShrinkByRadius ( const SizeType radius)

◆ ShrinkByRadius() [3/3]

template<unsigned int VImageDimension>
bool itk::ImageRegion< VImageDimension >::ShrinkByRadius ( OffsetValueType  radius)

Shrink an image region by the specified radius. The region can be shrunk uniformly in all dimension or can be shrink by different amounts in each direction. If the shrink operation fails because the radius is too large, this method returns false.

◆ Slice()

template<unsigned int VImageDimension>
SliceRegion itk::ImageRegion< VImageDimension >::Slice ( const unsigned int  dim) const

Slice a region, producing a region that is one dimension lower than the current region. Parameter "dim" specifies which dimension to remove.

Friends And Related Function Documentation

◆ ImageBase< VImageDimension >

template<unsigned int VImageDimension>
friend class ImageBase< VImageDimension >
friend

Friends of ImageRegion

Definition at line 419 of file itkImageRegion.h.

Member Data Documentation

◆ ImageDimension

template<unsigned int VImageDimension>
constexpr unsigned int itk::ImageRegion< VImageDimension >::ImageDimension = VImageDimension
staticconstexpr

Dimension of the image available at compile time.

Definition at line 102 of file itkImageRegion.h.

◆ m_Index

template<unsigned int VImageDimension>
IndexType itk::ImageRegion< VImageDimension >::m_Index = { { 0 } }
private

Definition at line 415 of file itkImageRegion.h.

Referenced by itk::ImageRegion< VDimension >::IsInside().

◆ m_Size

template<unsigned int VImageDimension>
SizeType itk::ImageRegion< VImageDimension >::m_Size = { { 0 } }
private

Definition at line 416 of file itkImageRegion.h.

Referenced by itk::ImageRegion< VDimension >::IsInside().

◆ SliceDimension

template<unsigned int VImageDimension>
constexpr unsigned int itk::ImageRegion< VImageDimension >::SliceDimension = ImageDimension - (ImageDimension > 1)
staticconstexpr

Dimension one lower than the image unless the image is one dimensional in which case the SliceDimension is also one dimensional.

Definition at line 106 of file itkImageRegion.h.


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