ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkRLEImage.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkRLEImage_h
19 #define itkRLEImage_h
20 
21 #include <itkImage.h>
22 #include <itkImageBase.h>
23 #include <utility> // std::pair
24 #include <vector>
25 
26 namespace itk
27 {
51 template< typename TPixel, unsigned int VImageDimension = 3, typename CounterType = unsigned short >
52 class RLEImage:
53  public itk::ImageBase< VImageDimension >
54 {
55 public:
56  ITK_DISALLOW_COPY_AND_ASSIGN(RLEImage);
57 
59  using Self = RLEImage;
64 
66  itkNewMacro( Self );
67 
69  itkTypeMacro( RLEImage, ImageBase );
70 
73  using PixelType = TPixel;
74 
75  using RLCounterType = CounterType;
76 
78  using ValueType = TPixel;
79 
82  using RLSegment = std::pair< CounterType, PixelType >;
83 
85  using RLLine = std::vector< RLSegment >;
86 
92 
93  // using IOPixelType = PixelType;
94 
99  static constexpr unsigned int ImageDimension = VImageDimension;
100 
104 
107 
109  using SizeType = typename Superclass::SizeType;
111 
114 
118 
123 
127 
130 
134  void
135  Allocate( bool initialize = false ) override;
136 
139  void
140  Initialize() override
141  {
142  // Call the superclass which should initialize the BufferedRegion ivar.
144  m_OnTheFlyCleanup = true;
145  m_Buffer = BufferType::New();
146  }
148 
151  void
152  FillBuffer( const TPixel& value );
153 
154  void
155  SetLargestPossibleRegion( const RegionType& region ) override
156  {
158  m_Buffer->SetLargestPossibleRegion( region.Slice( 0 ) );
159  }
160 
161  void
162  SetBufferedRegion( const RegionType& region ) override
163  {
165  m_Buffer->SetBufferedRegion( region.Slice( 0 ) );
166  }
167 
169 
170  void
171  SetRequestedRegion( const RegionType& region ) override
172  {
174  m_Buffer->SetRequestedRegion( region.Slice( 0 ) );
175  }
176 
182  void
183  SetPixel( const IndexType& index, const TPixel& value );
184 
189  int
190  SetPixel( RLLine& line, IndexValueType& segmentRemainder, SizeValueType& m_RealIndex, const TPixel& value );
191 
193  const TPixel &
194  GetPixel( const IndexType& index ) const;
195 
197  // TPixel & GetPixel(const IndexType & index);
198 
200  // TPixel & operator[](const IndexType & index)
201  // {
202  // return this->GetPixel(index);
203  // }
204 
207  const TPixel &
208  operator[]( const IndexType& index ) const
209  {
210  return this->GetPixel( index );
211  }
212 
213  unsigned int
215  {
216  // use the GetLength() method which works with variable length arrays,
217  // to make it work with as much pixel types as possible
218  PixelType p;
219 
220  return itk::NumericTraits< PixelType >::GetLength( p );
221  }
222 
224  using BufferType = typename itk::Image< RLLine, VImageDimension - 1 >;
225 
227  typename BufferType::Pointer
229  {
230  return m_Buffer;
231  }
232 
234  typename BufferType::Pointer
235  GetBuffer() const
236  {
237  return m_Buffer;
238  }
239 
241  static inline typename BufferType::IndexType
242  truncateIndex( const IndexType& index );
243 
246  void
247  CleanUp() const;
248 
251  bool
253  {
254  return m_OnTheFlyCleanup;
255  }
256 
259  void
260  SetOnTheFlyCleanup( bool value )
261  {
262  if ( value == m_OnTheFlyCleanup )
263  {
264  return;
265  }
266  m_OnTheFlyCleanup = value;
267  if ( m_OnTheFlyCleanup )
268  {
269  CleanUp(); // put the image into a clean state
270  }
271  }
273 
274 protected:
276  : itk::ImageBase< VImageDimension >(),
277  m_OnTheFlyCleanup( true )
278  {
279  m_Buffer = BufferType::New();
280  }
281 
282  void
283  PrintSelf( std::ostream& os, itk::Indent indent ) const override;
284 
285  ~RLEImage() override {}
291  void
293  {
295  }
296 
298  void
299  CleanUpLine( RLLine& line ) const;
300 
301 private:
302  bool m_OnTheFlyCleanup; // should same-valued segments be merged on the fly
303 
305  mutable typename BufferType::Pointer m_Buffer;
306 };
307 } // namespace itk
308 
309 #ifndef ITK_MANUAL_INSTANTIATION
310 #include "itkRLEImage.hxx"
311 #endif
312 
313 #endif // itkRLEImage_h
void ComputeIndexToPhysicalPointMatrices() override
Definition: itkRLEImage.h:292
void SetLargestPossibleRegion(const RegionType &region) override
Definition: itkRLEImage.h:155
BufferType::Pointer GetBuffer() const
Definition: itkRLEImage.h:235
static BufferType::IndexType truncateIndex(const IndexType &index)
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:151
BufferType::Pointer m_Buffer
Definition: itkRLEImage.h:305
TPixel PixelType
Definition: itkRLEImage.h:73
void PrintSelf(std::ostream &os, itk::Indent indent) const override
unsigned long SizeValueType
Definition: itkIntTypes.h:83
RLLine InternalPixelType
Definition: itkRLEImage.h:91
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Definition: itkImageBase.h:162
typename itk::Image< RLLine, VImageDimension-1 > BufferType
Definition: itkRLEImage.h:224
~RLEImage() override
Definition: itkRLEImage.h:285
bool GetOnTheFlyCleanup() const
Definition: itkRLEImage.h:252
ImageRegion< VImageDimension > RegionType
Definition: itkImageBase.h:145
An image region represents a structured region of data.
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
const TPixel & GetPixel(const IndexType &index) const
Get a pixel. SLOW! Better use iterators for pixel access.
Implements a weak reference to an object.
BufferType::Pointer GetBuffer()
Definition: itkRLEImage.h:228
void Initialize() override
Definition: itkRLEImage.h:140
std::vector< RLSegment > RLLine
Definition: itkRLEImage.h:85
const TPixel & operator[](const IndexType &index) const
Access a pixel. Chaning it changes the whole RLE segment!
Definition: itkRLEImage.h:208
Vector< SpacingValueType, VImageDimension > SpacingType
Definition: itkImageBase.h:152
signed long IndexValueType
Definition: itkIntTypes.h:90
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:68
Represent a n-dimensional offset between two n-dimensional indexes of n-dimensional image...
Definition: itkOffset.h:67
unsigned int GetNumberOfComponentsPerPixel() const override
Definition: itkRLEImage.h:214
void CleanUpLine(RLLine &line) const
typename SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:142
bool m_OnTheFlyCleanup
Definition: itkRLEImage.h:302
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:138
void FillBuffer(const TPixel &value)
void SetBufferedRegion(const RegionType &region) override
Definition: itkRLEImage.h:162
void SetOnTheFlyCleanup(bool value)
Definition: itkRLEImage.h:260
typename Superclass::IndexType IndexType
Definition: itkRLEImage.h:102
static constexpr unsigned int ImageDimension
Definition: itkRLEImage.h:99
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:137
Index< VImageDimension > IndexType
Definition: itkImageBase.h:132
Point< PointValueType, VImageDimension > PointType
Definition: itkImageBase.h:157
virtual void SetLargestPossibleRegion(const RegionType &region)
Base class for templated image classes.
Definition: itkImageBase.h:105
void Initialize() override
void Allocate(bool initialize=false) override
virtual void ComputeIndexToPhysicalPointMatrices()
CounterType RLCounterType
Definition: itkRLEImage.h:75
SliceRegion Slice(const unsigned int dim) const
TPixel ValueType
Definition: itkRLEImage.h:78
Control indentation during Print() invocation.
Definition: itkIndent.h:49
std::pair< CounterType, PixelType > RLSegment
Definition: itkRLEImage.h:82
Size< VImageDimension > SizeType
Definition: itkImageBase.h:141
Run-Length Encoded image. It saves memory for label images at the expense of processing times...
Definition: itkRLEImage.h:52
Base class for most ITK classes.
Definition: itkObject.h:60
typename IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:133
void SetRequestedRegion(const RegionType &region) override
Definition: itkRLEImage.h:171
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
void CleanUp() const
virtual void SetRequestedRegion(const RegionType &region)
virtual void SetBufferedRegion(const RegionType &region)