File: itkFlipImageFilter.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 (148 lines) | stat: -rw-r--r-- 5,504 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
/*=========================================================================
 *
 *  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 itkFlipImageFilter_h
#define itkFlipImageFilter_h

#include "itkImageToImageFilter.h"
#include "itkFixedArray.h"

namespace itk
{
/** \class FlipImageFilter
 * \brief Flips an image across user specified axes.
 *
 * FlipImageFilter flips an image across user specified axes.
 * The flip axes are set via method SetFlipAxes( array ) where
 * the input is a FixedArray<bool,ImageDimension>. The image
 * is flipped across axes for which array[i] is true.
 *
 * In terms of grid coordinates the image is flipped within
 * the LargestPossibleRegion of the input image. As such,
 * the LargestPossibleRegion of the output image is the same
 * as the input.
 *
 * In terms of geometric coordinates, the output origin
 * is such that the image is flipped with respect to the
 * coordinate axes.
 *
 * \ingroup GeometricTransform
 * \ingroup MultiThreaded
 * \ingroup Streamed
 * \ingroup ITKImageGrid
 *
 * \sphinx
 * \sphinxexample{Filtering/ImageGrid/FlipAnImageOverSpecifiedAxes,Flip An Image Over Specified Axes}
 * \endsphinx
 */
template <typename TImage>
class ITK_TEMPLATE_EXPORT FlipImageFilter : public ImageToImageFilter<TImage, TImage>
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(FlipImageFilter);

  /** Standard class type aliases. */
  using Self = FlipImageFilter;
  using Superclass = ImageToImageFilter<TImage, TImage>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;

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

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

  /** ImageDimension enumeration */
  static constexpr unsigned int ImageDimension = TImage::ImageDimension;

  /** Inherited types */
  using typename Superclass::InputImagePointer;
  using typename Superclass::InputImageConstPointer;
  using typename Superclass::OutputImagePointer;
  using typename Superclass::OutputImageRegionType;

  /** Index related types */
  using IndexType = typename TImage::IndexType;
  using IndexValueType = typename IndexType::IndexValueType;

  /** FlipAxesArray type */
  using FlipAxesArrayType = FixedArray<bool, Self::ImageDimension>;

  /** Set/Get the axis to be flipped. The image is flipped along axes
   * for which array[i] is true. Default is false. */
  itkSetMacro(FlipAxes, FlipAxesArrayType);
  itkGetConstMacro(FlipAxes, FlipAxesArrayType);

  /** Controls how the output origin is computed. If FlipAboutOrigin is
   * "On", the flip will occur about the origin of the axis, otherwise,
   * the flip will occur about the center of the axis. Default is "On". */
  itkBooleanMacro(FlipAboutOrigin);
  itkGetConstMacro(FlipAboutOrigin, bool);
  itkSetMacro(FlipAboutOrigin, bool);

  /** FlipImageFilter produces an image with different origin and
   * direction than the input image. As such, FlipImageFilter needs to
   * provide an implementation for GenerateOutputInformation() in
   * order to inform the pipeline execution model.
   * The output image meta information is obtained by permuting the input
   * image meta information. The original documentation of this method is
   * below.
   * \sa ProcessObject::GenerateOutputInformaton() */
  void
  GenerateOutputInformation() override;

  /** FlipImageFilter needs different input requested region than the output
   * requested region.  As such, FlipImageFilter needs to provide an
   * implementation for GenerateInputRequestedRegion() in order to inform the
   * pipeline execution model.
   * The required input requested region is obtained by permuting the index and
   * size of the output requested region.
   * \sa ProcessObject::GenerateInputRequestedRegion() */
  void
  GenerateInputRequestedRegion() override;

protected:
  FlipImageFilter();
  ~FlipImageFilter() override = default;
  void
  PrintSelf(std::ostream & os, Indent indent) const override;

  /** FlipImageFilter can be implemented as a multithreaded filter.
   * Therefore, this implementation provides a DynamicThreadedGenerateData() routine
   * which is called for each processing thread. The output image data is
   * allocated automatically by the superclass prior to calling
   * DynamicThreadedGenerateData().  DynamicThreadedGenerateData can only write to the
   * portion of the output image specified by the parameter
   * "outputRegionForThread"
   *
   * \sa ImageToImageFilter::ThreadedGenerateData(),
   *     ImageToImageFilter::GenerateData()  */
  void
  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;

private:
  FlipAxesArrayType m_FlipAxes{};
  bool              m_FlipAboutOrigin{ true };
};
} // end namespace itk

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

#endif