ITK  4.13.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:
57  typedef RLEImage Self;
62 
64  itkNewMacro( Self );
65 
67  itkTypeMacro( RLEImage, ImageBase );
68 
71  typedef TPixel PixelType;
72 
73  typedef CounterType RLCounterType;
74 
76  typedef TPixel ValueType;
77 
80  typedef std::pair< CounterType, PixelType > RLSegment;
81 
83  typedef std::vector< RLSegment > RLLine;
84 
90 
91  // typedef PixelType IOPixelType;
92 
97  itkStaticConstMacro( ImageDimension, unsigned int, VImageDimension );
98 
102 
105 
107  typedef typename Superclass::SizeType SizeType;
109 
112 
116 
121 
125 
128 
132  virtual void
133  Allocate( bool initialize = false );
134 
137  virtual void
139  {
140  // Call the superclass which should initialize the BufferedRegion ivar.
142  m_OnTheFlyCleanup = true;
144  }
146 
149  void
150  FillBuffer( const TPixel& value );
151 
152  virtual void
154  {
156  m_Buffer->SetLargestPossibleRegion( truncateRegion( region ) );
157  }
158 
159  virtual void
160  SetBufferedRegion( const RegionType& region )
161  {
163  m_Buffer->SetBufferedRegion( truncateRegion( region ) );
164  }
165 
167 
168  virtual void
170  {
172  m_Buffer->SetRequestedRegion( truncateRegion( region ) );
173  }
174 
180  void
181  SetPixel( const IndexType& index, const TPixel& value );
182 
187  int
188  SetPixel( RLLine& line, IndexValueType& segmentRemainder, SizeValueType& m_RealIndex, const TPixel& value );
189 
191  const TPixel &
192  GetPixel( const IndexType& index ) const;
193 
195  // TPixel & GetPixel(const IndexType & index);
196 
198  // TPixel & operator[](const IndexType & index)
199  // {
200  // return this->GetPixel(index);
201  // }
202 
205  const TPixel &
206  operator[]( const IndexType& index ) const
207  {
208  return this->GetPixel( index );
209  }
210 
211  virtual unsigned int
213  {
214  // use the GetLength() method which works with variable length arrays,
215  // to make it work with as much pixel types as possible
216  PixelType p;
217 
219  }
220 
222  typedef typename itk::Image< RLLine, VImageDimension - 1 > BufferType;
223 
225  typename BufferType::Pointer
227  {
228  return m_Buffer;
229  }
230 
232  typename BufferType::Pointer
233  GetBuffer() const
234  {
235  return m_Buffer;
236  }
237 
239  static inline typename BufferType::IndexType
240  truncateIndex( const IndexType& index );
241 
243  static inline typename BufferType::SizeType
244  truncateSize( const SizeType& size );
245 
247  static typename BufferType::RegionType
248  truncateRegion( const RegionType& region );
249 
252  void
253  CleanUp() const;
254 
257  bool
259  {
260  return m_OnTheFlyCleanup;
261  }
262 
265  void
266  SetOnTheFlyCleanup( bool value )
267  {
268  if ( value == m_OnTheFlyCleanup )
269  {
270  return;
271  }
272  m_OnTheFlyCleanup = value;
273  if ( m_OnTheFlyCleanup )
274  {
275  CleanUp(); // put the image into a clean state
276  }
277  }
279 
280 protected:
282  : itk::ImageBase< VImageDimension >(),
283  m_OnTheFlyCleanup( true )
284  {
285  // m_OnTheFlyCleanup = true;
287  }
288 
289  void
290  PrintSelf( std::ostream& os, itk::Indent indent ) const;
291 
292  virtual ~RLEImage() {}
298  virtual void
300  {
302  }
303 
305  void
306  CleanUpLine( RLLine& line ) const;
307 
308 private:
309  bool m_OnTheFlyCleanup; // should same-valued segments be merged on the fly
310 
311  RLEImage( const Self & ); // purposely not implemented
312  void
313  operator=( const Self & ); // purposely not implemented
314 
316  mutable typename BufferType::Pointer m_Buffer;
317 };
318 } // namespace itk
319 
320 #ifndef ITK_MANUAL_INSTANTIATION
321 #include "itkRLEImage.hxx"
322 #endif
323 
324 #endif // itkRLEImage_h
Superclass::PointType PointType
Definition: itkRLEImage.h:124
Superclass::RegionType RegionType
Definition: itkImage.h:137
itk::Image< RLLine, VImageDimension-1 > BufferType
Definition: itkRLEImage.h:222
RLLine InternalPixelType
Definition: itkRLEImage.h:89
Index< VImageDimension > IndexType
Definition: itkImageBase.h:139
virtual void SetRequestedRegion(const RegionType &region)
Definition: itkRLEImage.h:169
BufferType::Pointer GetBuffer() const
Definition: itkRLEImage.h:233
SmartPointer< Self > Pointer
Definition: itkImage.h:81
Superclass::OffsetValueType OffsetValueType
Definition: itkRLEImage.h:127
static BufferType::IndexType truncateIndex(const IndexType &index)
Superclass::SpacingType SpacingType
Definition: itkRLEImage.h:119
BufferType::Pointer m_Buffer
Definition: itkRLEImage.h:316
Superclass::OffsetType OffsetType
Definition: itkRLEImage.h:104
virtual void ComputeIndexToPhysicalPointMatrices()
Definition: itkRLEImage.h:299
itk::WeakPointer< const Self > ConstWeakPointer
Definition: itkRLEImage.h:61
signed long IndexValueType
Definition: itkIntTypes.h:150
std::pair< CounterType, PixelType > RLSegment
Definition: itkRLEImage.h:80
Superclass::SpacingValueType SpacingValueType
Definition: itkRLEImage.h:120
bool GetOnTheFlyCleanup() const
Definition: itkRLEImage.h:258
An image region represents a structured region of data.
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:158
Superclass::SizeType SizeType
Definition: itkRLEImage.h:107
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
TPixel PixelType
Definition: itkRLEImage.h:67
static Pointer New()
Superclass::DirectionType DirectionType
Definition: itkRLEImage.h:111
const TPixel & GetPixel(const IndexType &index) const
Get a pixel. SLOW! Better use iterators for pixel access.
Implements a weak reference to an object.
virtual void Initialize() override
unsigned long SizeValueType
Definition: itkIntTypes.h:143
static const unsigned int ImageDimension
Definition: itkRLEImage.h:97
Size< VImageDimension > SizeType
Definition: itkImageBase.h:148
Point< PointValueType, VImageDimension > PointType
Definition: itkImageBase.h:164
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Definition: itkImageBase.h:169
BufferType::Pointer GetBuffer()
Definition: itkRLEImage.h:226
virtual void Allocate(bool initialize=false)
Superclass::RegionType RegionType
Definition: itkRLEImage.h:115
const TPixel & operator[](const IndexType &index) const
Access a pixel. Chaning it changes the whole RLE segment!
Definition: itkRLEImage.h:206
static BufferType::SizeType truncateSize(const SizeType &size)
Superclass::IndexType IndexType
Definition: itkRLEImage.h:100
static BufferType::RegionType truncateRegion(const RegionType &region)
static unsigned int GetLength()
RLEImage Self
Definition: itkRLEImage.h:57
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:144
virtual ~RLEImage()
Definition: itkRLEImage.h:292
Superclass::IndexType IndexType
Definition: itkImage.h:119
Vector< SpacingValueType, VImageDimension > SpacingType
Definition: itkImageBase.h:159
Superclass::SizeType SizeType
Definition: itkImage.h:126
void CleanUpLine(RLLine &line) const
bool m_OnTheFlyCleanup
Definition: itkRLEImage.h:309
void FillBuffer(const TPixel &value)
void SetOnTheFlyCleanup(bool value)
Definition: itkRLEImage.h:266
Superclass::SizeValueType SizeValueType
Definition: itkRLEImage.h:108
itk::SmartPointer< const Self > ConstPointer
Definition: itkRLEImage.h:60
SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:149
OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:145
virtual void SetLargestPossibleRegion(const RegionType &region)
Definition: itkRLEImage.h:153
virtual unsigned int GetNumberOfComponentsPerPixel() const
Definition: itkRLEImage.h:212
virtual void SetLargestPossibleRegion(const RegionType &region)
Base class for templated image classes.
Definition: itkImageBase.h:114
virtual void ComputeIndexToPhysicalPointMatrices()
void operator=(const Self &)
itk::SmartPointer< Self > Pointer
Definition: itkRLEImage.h:59
Control indentation during Print() invocation.
Definition: itkIndent.h:49
TPixel ValueType
Definition: itkRLEImage.h:76
void PrintSelf(std::ostream &os, itk::Indent indent) const
itk::ImageBase< VImageDimension > Superclass
Definition: itkRLEImage.h:58
virtual void SetBufferedRegion(const RegionType &region)
Definition: itkRLEImage.h:160
Run-Length Encoded image. It saves memory for label images at the expense of processing times...
Definition: itkRLEImage.h:52
virtual void Initialize()
Definition: itkRLEImage.h:138
std::vector< RLSegment > RLLine
Definition: itkRLEImage.h:83
Base class for all data objects in ITK.
IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:140
Superclass::IndexValueType IndexValueType
Definition: itkRLEImage.h:101
Templated n-dimensional image class.
Definition: itkImage.h:75
void CleanUp() const
virtual void SetRequestedRegion(const RegionType &region)
CounterType RLCounterType
Definition: itkRLEImage.h:73
virtual void SetBufferedRegion(const RegionType &region)