File: itkImage.h

package info (click to toggle)
insighttoolkit5 5.4.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 704,384 kB
  • sloc: cpp: 783,592; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 464; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (404 lines) | stat: -rw-r--r-- 13,779 bytes parent folder | download
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
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
/*=========================================================================
 *
 *  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 itkImage_h
#define itkImage_h

#include "itkImageRegion.h"
#include "itkImportImageContainer.h"
#include "itkDefaultPixelAccessor.h"
#include "itkDefaultPixelAccessorFunctor.h"
#include "itkPoint.h"
#include "itkFixedArray.h"
#include "itkWeakPointer.h"
#include "itkNeighborhoodAccessorFunctor.h"

#include <type_traits> // For is_same

namespace itk
{
/** \class Image
 *  \brief Templated n-dimensional image class.
 *
 * Images are templated over a pixel type (modeling the dependent
 * variables), and a dimension (number of independent variables).  The
 * container for the pixel data is the ImportImageContainer.
 *
 * Within the pixel container, images are modeled as arrays, defined by a
 * start index and a size.
 *
 * The superclass of Image, ImageBase, defines the geometry of the
 * image in terms of where the image sits in physical space, how the
 * image is oriented in physical space, the size of a pixel, and the
 * extent of the image itself.  ImageBase provides the methods to
 * convert between the index and physical space coordinate frames.
 *
 * Pixels can be accessed directly using the SetPixel() and GetPixel()
 * methods or can be accessed via iterators that define the region of
 * the image they traverse.
 *
 * The pixel type may be one of the native types; a Insight-defined
 * class type such as Vector; or a user-defined type. Note that
 * depending on the type of pixel that you use, the process objects
 * (i.e., those filters processing data objects) may not operate on
 * the image and/or pixel type. This becomes apparent at compile-time
 * because operator overloading (for the pixel type) is not supported.
 *
 * The data in an image is arranged in a 1D array as [][][][slice][row][col]
 * with the column index varying most rapidly.  The Index type reverses
 * the order so that with Index[0] = col, Index[1] = row, Index[2] = slice,
 * ...
 *
 * \sa ImageBase
 * \sa ImageContainerInterface
 *
 * \ingroup ImageObjects
 * \ingroup ITKCommon
 *
 * \sphinx
 * \sphinxexample{Core/Common/SetPixelValueInOneImage,Set Pixel Value In One Image}
 * \sphinxexample{Core/Common/GetImageSize,Get Image Size}
 * \sphinxexample{Core/Common/SortITKIndex,Sort ITK Index}
 * \sphinxexample{Core/Common/ReturnObjectFromFunction,Return Object From Function}
 * \sphinxexample{Core/Common/CreateAnotherInstanceOfAnImage,Create Another Instance Of An Image}
 * \sphinxexample{Core/Common/PassImageToFunction,Pass Image To Function}
 * \sphinxexample{Core/Common/DeepCopyImage,Deep Copy Image}
 * \sphinxexample{Core/Common/ThrowException,Throw Exception}
 * \sphinxexample{Core/Common/GetOrSetMemberVariableOfITKClass,Get Or Set Member Variable Of ITK Class}
 * \sphinxexample{Core/Common/MiniPipeline,Mini Pipeline}
 * \sphinxexample{Core/Common/CheckIfModuleIsPresent,Check If Module Is Present}
 * \sphinxexample{Core/Common/DisplayImage,Display Image}
 * \endsphinx
 */
template <typename TPixel, unsigned int VImageDimension = 2>
class ITK_TEMPLATE_EXPORT Image : public ImageBase<VImageDimension>
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(Image);

  /** Standard class type aliases */
  using Self = Image;
  using Superclass = ImageBase<VImageDimension>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;
  using ConstWeakPointer = WeakPointer<const Self>;

  /** Method for creation through the object factory. */
  itkNewMacro(Self);

  /** \see LightObject::GetNameOfClass() */
  itkOverrideGetNameOfClassMacro(Image);

  /** Pixel type alias support. Used to declare pixel type in filters
   * or other operations. */
  using PixelType = TPixel;

  /** Typedef alias for PixelType */
  using ValueType = TPixel;

  /** Internal Pixel representation. Used to maintain a uniform API
   * with Image Adaptors and allow to keep a particular internal
   * representation of data while showing a different external
   * representation. */
  using InternalPixelType = TPixel;

  using IOPixelType = PixelType;

  /** Accessor type that convert data between internal and external
   *  representations.  */
  using AccessorType = DefaultPixelAccessor<PixelType>;
  using AccessorFunctorType = DefaultPixelAccessorFunctor<Self>;

  /** Typedef for the functor used to access a neighborhood of pixel
   * pointers. */
  using NeighborhoodAccessorFunctorType = NeighborhoodAccessorFunctor<Self>;

  /** Type of image dimension */
  using typename Superclass::ImageDimensionType;

  /** Index type alias support. An index is used to access pixel values. */
  using typename Superclass::IndexType;
  using typename Superclass::IndexValueType;

  /** Offset type alias support. An offset is used to access pixel values. */
  using typename Superclass::OffsetType;

  /** Size type alias support. A size is used to define region bounds. */
  using typename Superclass::SizeType;
  using typename Superclass::SizeValueType;

  /** Container used to store pixels in the image. */
  using PixelContainer = ImportImageContainer<SizeValueType, PixelType>;

  /** Direction type alias support. A matrix of direction cosines. */
  using typename Superclass::DirectionType;

  /** Region type alias support. A region is used to specify a subset of an image.
   */
  using typename Superclass::RegionType;

  /** Spacing type alias support.  Spacing holds the size of a pixel.  The
   * spacing is the geometric distance between image samples. */
  using typename Superclass::SpacingType;
  using typename Superclass::SpacingValueType;

  /** Origin type alias support.  The origin is the geometric coordinates
   * of the index (0,0). */
  using typename Superclass::PointType;

  /** A pointer to the pixel container. */
  using PixelContainerPointer = typename PixelContainer::Pointer;
  using PixelContainerConstPointer = typename PixelContainer::ConstPointer;

  /** Offset type alias (relative position between indices) */
  using typename Superclass::OffsetValueType;

  /**
   * example usage:
   * using OutputImageType = typename ImageType::template Rebind< float >::Type;
   *
   * \deprecated Use RebindImageType instead
   */
  template <typename UPixelType, unsigned int VUImageDimension = VImageDimension>
  struct Rebind
  {
    using Type = itk::Image<UPixelType, VUImageDimension>;
  };


  template <typename UPixelType, unsigned int VUImageDimension = VImageDimension>
  using RebindImageType = itk::Image<UPixelType, VUImageDimension>;


  /** Allocate the image memory. The size of the image must
   * already be set, e.g. by calling SetRegions(). */
  void
  Allocate(bool initializePixels = false) override;

  /** Restore the data object to its initial state. This means releasing
   * memory. */
  void
  Initialize() override;

  /** Fill the image buffer with a value.  Be sure to call Allocate()
   * first. */
  void
  FillBuffer(const TPixel & value);

  /** \brief Set a pixel value.
   *
   * Allocate() needs to have been called first -- for efficiency,
   * this function does not check that the image has actually been
   * allocated yet. */
  void
  SetPixel(const IndexType & index, const TPixel & value)
  {
    OffsetValueType offset = this->FastComputeOffset(index);
    (*m_Buffer)[offset] = value;
  }

  /** \brief Get a pixel (read only version).
   *
   * For efficiency, this function does not check that the
   * image has actually been allocated yet. */
  const TPixel &
  GetPixel(const IndexType & index) const
  {
    OffsetValueType offset = this->FastComputeOffset(index);
    return ((*m_Buffer)[offset]);
  }

  /** \brief Get a reference to a pixel (e.g. for editing).
   *
   * For efficiency, this function does not check that the
   * image has actually been allocated yet. */
  TPixel &
  GetPixel(const IndexType & index)
  {
    OffsetValueType offset = this->FastComputeOffset(index);
    return ((*m_Buffer)[offset]);
  }

  /** \brief Access a pixel. This version can be an lvalue.
   *
   * For efficiency, this function does not check that the
   * image has actually been allocated yet. */
  TPixel & operator[](const IndexType & index) { return this->GetPixel(index); }

  /** \brief Access a pixel. This version can only be an rvalue.
   *
   * For efficiency, this function does not check that the
   * image has actually been allocated yet. */
  const TPixel & operator[](const IndexType & index) const { return this->GetPixel(index); }

  /** Return a pointer to the beginning of the buffer.  This is used by
   * the image iterator class. */
  virtual TPixel *
  GetBufferPointer()
  {
    return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
  }
  virtual const TPixel *
  GetBufferPointer() const
  {
    return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
  }

  /** Return a pointer to the container. */
  PixelContainer *
  GetPixelContainer()
  {
    return m_Buffer.GetPointer();
  }

  const PixelContainer *
  GetPixelContainer() const
  {
    return m_Buffer.GetPointer();
  }

  /** Set the container to use. Note that this does not cause the
   * DataObject to be modified. */
  void
  SetPixelContainer(PixelContainer * container);

  /** Graft the data and information from one image to another. This
   * is a convenience method to setup a second image with all the meta
   * information of another image and use the same pixel
   * container. Note that this method is different than just using two
   * SmartPointers to the same image since separate DataObjects are
   * still maintained. This method is similar to
   * ImageSource::GraftOutput(). The implementation in ImageBase
   * simply calls CopyInformation() and copies the region ivars.
   * The implementation here refers to the superclass' implementation
   * and then copies over the pixel container. */
  virtual void
  Graft(const Self * image);

  /** Return the Pixel Accessor object */
  AccessorType
  GetPixelAccessor()
  {
    return AccessorType();
  }

  /** Return the Pixel Accesor object */
  const AccessorType
  GetPixelAccessor() const
  {
    return AccessorType();
  }

  /** Return the NeighborhoodAccessor functor */
  NeighborhoodAccessorFunctorType
  GetNeighborhoodAccessor()
  {
    return NeighborhoodAccessorFunctorType();
  }

  /** Return the NeighborhoodAccessor functor */
  const NeighborhoodAccessorFunctorType
  GetNeighborhoodAccessor() const
  {
    return NeighborhoodAccessorFunctorType();
  }

  unsigned int
  GetNumberOfComponentsPerPixel() const override;

  /** Returns (image1 == image2).
   * \note `operator==` and `operator!=` are defined as function templates
   * (rather than as non-templates), just to allow template instantiation of
   * `itk::Image` for non-EqualityComparable pixel types.
   */
  template <typename TEqualityComparable>
  friend std::enable_if_t<std::is_same_v<TEqualityComparable, TPixel>, bool>
  operator==(const Image<TEqualityComparable, VImageDimension> & lhs,
             const Image<TEqualityComparable, VImageDimension> & rhs)
  {
    if ((lhs.GetBufferedRegion() != rhs.GetBufferedRegion()) || (lhs.m_Spacing != rhs.m_Spacing) ||
        (lhs.m_Origin != rhs.m_Origin) || (lhs.m_Direction != rhs.m_Direction) ||
        (lhs.m_InverseDirection != rhs.m_InverseDirection))
    {
      return false;
    }

    if (lhs.m_Buffer == rhs.m_Buffer)
    {
      return true;
    }

    if ((lhs.m_Buffer == nullptr) || (rhs.m_Buffer == nullptr))
    {
      return false;
    }

    auto & lhsBuffer = *(lhs.m_Buffer);
    auto & rhsBuffer = *(rhs.m_Buffer);

    const auto bufferSize = lhsBuffer.Size();

    if (bufferSize != rhsBuffer.Size())
    {
      return false;
    }

    const TEqualityComparable * const lhsBufferPointer = lhsBuffer.GetBufferPointer();
    const TEqualityComparable * const rhsBufferPointer = rhsBuffer.GetBufferPointer();

    return ((lhsBufferPointer == rhsBufferPointer) ||
            std::equal(lhsBufferPointer, lhsBufferPointer + bufferSize, rhsBufferPointer));
  }

  /** Returns (image1 != image2). */
  template <typename TEqualityComparable>
  friend std::enable_if_t<std::is_same_v<TEqualityComparable, TPixel>, bool>
  operator!=(const Image<TEqualityComparable, VImageDimension> & lhs,
             const Image<TEqualityComparable, VImageDimension> & rhs)
  {
    return !(lhs == rhs);
  }

protected:
  Image() = default;
  void
  PrintSelf(std::ostream & os, Indent indent) const override;
  void
  Graft(const DataObject * data) override;

  ~Image() override = default;

  /** Compute helper matrices used to transform Index coordinates to
   * PhysicalPoint coordinates and back. This method is virtual and will be
   * overloaded in derived classes in order to provide backward compatibility
   * behavior in classes that did not used to take image orientation into
   * account.  */
  void
  ComputeIndexToPhysicalPointMatrices() override;
  using Superclass::Graft;

private:
  /** Memory for the current buffer. */
  PixelContainerPointer m_Buffer{ PixelContainer::New() };
};
} // end namespace itk

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

#endif