RTK  2.6.0
Reconstruction Toolkit
rtkConstantImageSource.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright RTK 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  * https://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 
19 #ifndef rtkConstantImageSource_h
20 #define rtkConstantImageSource_h
21 
22 #include "rtkConfiguration.h"
23 #include "rtkMacro.h"
24 
25 #include <itkImageSource.h>
26 #include <itkNumericTraits.h>
28 #include <itkVectorImage.h>
29 
30 namespace rtk
31 {
50 template <typename TOutputImage>
51 class ITK_TEMPLATE_EXPORT ConstantImageSource : public itk::ImageSource<TOutputImage>
52 {
53 public:
54  ITK_DISALLOW_COPY_AND_MOVE(ConstantImageSource);
55 
61 
63  using OutputImageType = TOutputImage;
64 
66  using OutputImagePixelType = typename TOutputImage::PixelType;
67 
69  using OutputImageRegionType = typename TOutputImage::RegionType;
70 
72  itkOverrideGetNameOfClassMacro(ConstantImageSource);
73 
75  itkNewMacro(Self);
76 
78  using SizeType = typename TOutputImage::SizeType;
79  using IndexType = typename TOutputImage::IndexType;
80  using SpacingType = typename TOutputImage::SpacingType;
81  using PointType = typename TOutputImage::PointType;
83  typedef SizeValueType SizeValueArrayType[TOutputImage::ImageDimension];
84  using SpacingValueType = typename TOutputImage::SpacingValueType;
85  typedef SpacingValueType SpacingValueArrayType[TOutputImage::ImageDimension];
86  using PointValueType = typename TOutputImage::PointValueType;
87  typedef PointValueType PointValueArrayType[TOutputImage::ImageDimension];
88  using DirectionType = typename TOutputImage::DirectionType;
89 
92  itkGetMacro(Size, SizeType);
93  virtual void
94  SetSize(SizeValueArrayType sizeArray);
95  virtual const SizeValueType *
96  GetSize() const;
98 
100  itkSetMacro(Spacing, SpacingType);
101  itkGetMacro(Spacing, SpacingType);
103 
105  itkSetMacro(Origin, PointType);
106  itkGetMacro(Origin, PointType);
108 
110  itkSetMacro(Direction, DirectionType);
111  itkGetMacro(Direction, DirectionType);
113 
116  itkGetMacro(Index, IndexType);
118 
121  itkGetConstMacro(Constant, OutputImagePixelType);
123 
125  void
126  SetInformationFromImage(const itk::ImageBase<TOutputImage::ImageDimension> * image);
127 
128 protected:
130  ~ConstantImageSource() override;
131  void
132  PrintSelf(std::ostream & os, itk::Indent indent) const override;
133 
134  void
135  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
136 
137  void
138  GenerateOutputInformation() override;
139 
145 
147 };
148 
149 } // end namespace rtk
150 
151 #ifndef ITK_MANUAL_INSTANTIATION
152 # include "rtkConstantImageSource.hxx"
153 #endif
154 
155 #endif
typename TOutputImage::SizeType SizeType
typename SizeType::SizeValueType SizeValueType
Generate an n-dimensional image with constant pixel values.
typename TOutputImage::SpacingType SpacingType
typename TOutputImage::DirectionType DirectionType
typename TOutputImage::PointValueType PointValueType
#define itkSetMacro(name, type)
typename TOutputImage::RegionType OutputImageRegionType
typename TOutputImage::IndexType IndexType
OutputImagePixelType m_Constant
typename TOutputImage::PointType PointType
typename TOutputImage::PixelType OutputImagePixelType
itk::SizeValueType SizeValueType
typename TOutputImage::SpacingValueType SpacingValueType