1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136
|
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#ifndef itkImageRegionIteratorWithIndex_h
#define itkImageRegionIteratorWithIndex_h
#include "itkImageRegionConstIteratorWithIndex.h"
#include "itkImageIteratorWithIndex.h"
namespace itk
{
/** \class ImageRegionIteratorWithIndex
* \brief A multi-dimensional iterator templated over image type that walks
* pixels within a region and is specialized to keep track of its image index
* location.
*
* This class is a specialization of ImageRegionConstIteratorWithIndex that
* adds write-access (the Set() method). Please see
* ImageRegionConstIteratorWithIndex for more information.
*
* \par MORE INFORMATION
* For a complete description of the ITK Image Iterators and their API, please
* see the Iterators chapter in the ITK Software Guide. The ITK Software Guide
* is available in print and as a free .pdf download from https://www.itk.org.
*
* \ingroup ImageIterators
*
* \sa ImageConstIterator \sa ConditionalConstIterator
* \sa ConstNeighborhoodIterator \sa ConstShapedNeighborhoodIterator
* \sa ConstSliceIterator \sa CorrespondenceDataStructureIterator
* \sa FloodFilledFunctionConditionalConstIterator
* \sa FloodFilledImageFunctionConditionalConstIterator
* \sa FloodFilledImageFunctionConditionalIterator
* \sa FloodFilledSpatialFunctionConditionalConstIterator
* \sa FloodFilledSpatialFunctionConditionalIterator
* \sa ImageConstIterator \sa ImageConstIteratorWithIndex
* \sa ImageIterator \sa ImageIteratorWithIndex
* \sa ImageLinearConstIteratorWithIndex \sa ImageLinearIteratorWithIndex
* \sa ImageRandomConstIteratorWithIndex \sa ImageRandomIteratorWithIndex
* \sa ImageRegionConstIterator \sa ImageRegionConstIteratorWithIndex
* \sa ImageRegionExclusionConstIteratorWithIndex
* \sa ImageRegionExclusionIteratorWithIndex
* \sa ImageRegionIterator \sa ImageRegionIteratorWithIndex
* \sa ImageRegionReverseConstIterator \sa ImageRegionReverseIterator
* \sa ImageReverseConstIterator \sa ImageReverseIterator
* \sa ImageSliceConstIteratorWithIndex \sa ImageSliceIteratorWithIndex
* \sa NeighborhoodIterator \sa PathConstIterator \sa PathIterator
* \sa ShapedNeighborhoodIterator \sa SliceIterator
* \sa ImageConstIteratorWithIndex
* \sa ImageRegionRange
* \sa ImageRegionIndexRange
* \ingroup ITKCommon
*
* \sphinx
* \sphinxexample{Core/Common/IterateRegionWithAccessToIndexWithWriteAccess,Iterate Region In Image With Access To
* Current Index With Write Access} \endsphinx
*/
template <typename TImage>
class ITK_TEMPLATE_EXPORT ImageRegionIteratorWithIndex : public ImageRegionConstIteratorWithIndex<TImage>
{
public:
/** Standard class type aliases. */
using Self = ImageRegionIteratorWithIndex;
using Superclass = ImageRegionConstIteratorWithIndex<TImage>;
/** Types inherited from the Superclass */
using typename Superclass::IndexType;
using typename Superclass::SizeType;
using typename Superclass::OffsetType;
using typename Superclass::RegionType;
using typename Superclass::ImageType;
using typename Superclass::PixelContainer;
using typename Superclass::PixelContainerPointer;
using typename Superclass::InternalPixelType;
using typename Superclass::PixelType;
using typename Superclass::AccessorType;
/** Default constructor. Needed since we provide a cast constructor. */
ImageRegionIteratorWithIndex() = default;
/** Constructor establishes an iterator to walk a particular image and a
* particular region of that image. */
ImageRegionIteratorWithIndex(TImage * ptr, const RegionType & region);
/** Constructor that can be used to cast from an ImageIterator to an
* ImageRegionIteratorWithIndex. Many routines return an ImageIterator, but for a
* particular task, you may want an ImageRegionIteratorWithIndex. Rather than
* provide overloaded APIs that return different types of Iterators, itk
* returns ImageIterators and uses constructors to cast from an
* ImageIterator to a ImageRegionIteratorWithIndex. */
ImageRegionIteratorWithIndex(const ImageIteratorWithIndex<TImage> & it);
/** Set the pixel value */
void
Set(const PixelType & value) const
{
this->m_PixelAccessorFunctor.Set(*(const_cast<InternalPixelType *>(this->m_Position)), value);
}
/** Return a reference to the pixel.
* This method will provide the fastest access to pixel
* data, but it will NOT support ImageAdaptors. */
PixelType &
Value()
{
return *(const_cast<InternalPixelType *>(this->m_Position));
}
protected:
/** The construction from a const iterator is declared protected
in order to enforce const correctness. */
ImageRegionIteratorWithIndex(const ImageRegionConstIteratorWithIndex<TImage> & it);
Self &
operator=(const ImageRegionConstIteratorWithIndex<TImage> & it);
};
} // end namespace itk
#ifndef ITK_MANUAL_INSTANTIATION
# include "itkImageRegionIteratorWithIndex.hxx"
#endif
#endif
|