File: itkLineConstIterator.h

package info (click to toggle)
insighttoolkit5 5.4.4-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 704,404 kB
  • sloc: cpp: 783,697; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 461; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (199 lines) | stat: -rw-r--r-- 6,205 bytes parent folder | download | duplicates (2)
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
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
/*=========================================================================
 *
 *  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 itkLineConstIterator_h
#define itkLineConstIterator_h

#include "itkIndex.h"
#include "itkImage.h"

namespace itk
{
/**
 * \class LineConstIterator
 * \brief An iterator that walks a Bresenham line through an ND image
 *        with read-only access to pixels.
 *
 * LineConstIterator is an iterator that walks a Bresenham line
 * through an image.  The iterator is constructed similar to other
 * image iterators, except instead of specifying a region to
 * traverse, you specify two indices. The interval specified by
 * the two indices is closed.  So, a line iterator specified with
 * the same start and end index will visit exactly one pixel.
 *
   \code
   LineConstIterator<ImageType> it(image, I1, I2);
   while (!it.IsAtEnd())
   {
      // visits at least 1 pixel
   }
   \endcode
 *
 * \author Benjamin King, Experimentelle Radiologie, Medizinische
 * Hochschule Hannover.
 *
 * \ingroup ITKCommon
 *
 * \sphinx
 * \sphinxexample{Core/Common/IterateLineThroughImageWithoutWriteAccess,Iterate Line Through Image Without Write Access}
 * \endsphinx
 */
template <typename TImage>
class ITK_TEMPLATE_EXPORT LineConstIterator
{
public:
  /** Standard class type aliases. */
  using Self = LineConstIterator;

  /** Dimension of the image that the iterator walks.  This constant is needed so
   * that functions that are templated over image iterator type (as opposed to
   * being templated over pixel type and dimension) can have compile time
   * access to the dimension of the image that the iterator walks. */
  static constexpr unsigned int ImageIteratorDimension = TImage::ImageDimension;

  /** Index type alias support */
  using IndexType = typename TImage::IndexType;

  /** Offset type alias support */
  using OffsetType = typename TImage::OffsetType;

  /** Size type alias support */
  using SizeType = typename TImage::SizeType;

  /** Region type alias support */
  using RegionType = typename TImage::RegionType;

  /** Spacing type alias support */
  using SpacingType = typename TImage::SpacingType;

  /** Origin type alias support */
  using PointType = typename TImage::PointType;

  /** Image type alias support */
  using ImageType = TImage;

  /** PixelContainer type alias support Used to refer to the container for
   * the pixel data. While this was already typedef'ed in the superclass,
   * it needs to be redone here for this subclass to compile properly with gcc. */
  using PixelContainer = typename TImage::PixelContainer;
  using PixelContainerPointer = typename PixelContainer::Pointer;

  /** Internal Pixel Type */
  using InternalPixelType = typename TImage::InternalPixelType;

  /** External Pixel Type */
  using PixelType = typename TImage::PixelType;

  /**  Accessor type that convert data between internal and external
   *  representations. */
  using AccessorType = typename TImage::AccessorType;

  /** \see LightObject::GetNameOfClass() */
  itkVirtualGetNameOfClassMacro(LineConstIterator);

  /** Get the dimension (size) of the index. */
  static unsigned int
  GetImageIteratorDimension()
  {
    return TImage::ImageDimension;
  }

  /** Get the index. This provides a read only reference to the index. */
  const IndexType
  GetIndex()
  {
    return m_CurrentImageIndex;
  }

  /** Get the pixel value */
  const PixelType
  Get() const
  {
    return m_Image->GetPixel(m_CurrentImageIndex);
  }

  /** Is the iterator at the end of the line? */
  bool
  IsAtEnd() const
  {
    return m_IsAtEnd;
  }

  /** Move an iterator to the beginning of the line. */
  void
  GoToBegin();

  /** Walk forward along the line to the next index in the image. */
  void
  operator++();

  /** operator= is provided to make sure the handle to the image is properly
   * reference counted. */
  Self &
  operator=(const Self & it);

  /** Constructor establishes an iterator to walk along a line */
  LineConstIterator(const ImageType * imagePtr, const IndexType & firstIndex, const IndexType & lastIndex);

  /** Default Destructor. */
  virtual ~LineConstIterator() = default;

protected: // made protected so other iterators can access
  /** Smart pointer to the source image. */
  typename ImageType::ConstWeakPointer m_Image{};

  /** Region type to iterate over. */
  RegionType m_Region{};

  /** Is the iterator at the end of its walk? */
  bool m_IsAtEnd{};

  /** Start, end and current ND index position in the image of the line */
  IndexType m_CurrentImageIndex{};
  IndexType m_StartIndex{};
  IndexType m_LastIndex{};
  IndexType m_EndIndex{}; // one past the end of the line in the m_MainDirection

  /** Variables that drive the Bresenham-Algorithm */
  // The dimension with the largest difference between start and end
  unsigned int m_MainDirection{};

  // Accumulated error for the other dimensions
  IndexType m_AccumulateError{};

  // Increment for the error for each step. Two times the difference between
  // start and end
  IndexType m_IncrementError{};

  // If enough is accumulated for a dimension, the index has to be
  // incremented. Will be the number of pixels in the line
  IndexType m_MaximalError{};

  // Direction of increment. -1 or 1
  IndexType m_OverflowIncrement{};

  // After an overflow, the accumulated error is reduced again. Will be
  // two times the number of pixels in the line
  IndexType m_ReduceErrorAfterIncrement{};
};
} // end namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#  include "itkLineConstIterator.hxx"
#endif

#endif