RTK  2.6.0
Reconstruction Toolkit
rtkAdditiveGaussianNoiseImageFilter.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 /*=========================================================================
20  *
21  * Copyright NumFOCUS
22  *
23  * Licensed under the Apache License, Version 2.0 (the "License");
24  * you may not use this file except in compliance with the License.
25  * You may obtain a copy of the License at
26  *
27  * https://www.apache.org/licenses/LICENSE-2.0.txt
28  *
29  * Unless required by applicable law or agreed to in writing, software
30  * distributed under the License is distributed on an "AS IS" BASIS,
31  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
32  * See the License for the specific language governing permissions and
33  * limitations under the License.
34  *
35  *=========================================================================*/
36 #ifndef rtkAdditiveGaussianNoiseImageFilter_h
37 #define rtkAdditiveGaussianNoiseImageFilter_h
38 
39 #include <itkImageToImageFilter.h>
42 
43 #include "rtkMacro.h"
44 
45 namespace rtk
46 {
47 
56 template <class TPixel>
57 class ITK_TEMPLATE_EXPORT NormalVariateNoiseFunctor
58 {
59 public:
61  {
62  m_Mean = 0.0;
63  m_StandardDeviation = 1.0;
64 
66 
67  this->SetSeed(42);
68  }
69 
70  float
71  GetMean() const
72  {
73  return m_Mean;
74  }
75 
76  void
77  SetMean(float mean)
78  {
79  m_Mean = mean;
80  }
81 
82  float
84  {
85  return m_StandardDeviation;
86  }
87 
88  void
89  SetStandardDeviation(float stddev)
90  {
91  m_StandardDeviation = stddev;
92  }
93 
94  void
95  SetSeed(unsigned long seed)
96  {
97  m_Generator->Initialize(seed);
98  }
99 
100  void
101  SetOutputMinimum(TPixel min)
102  {
103  m_OutputMinimum = min;
104  }
105 
106  void
107  SetOutputMaximum(TPixel max)
108  {
109  m_OutputMaximum = max;
110  }
111 
112  TPixel
114  {
115  return m_OutputMinimum;
116  }
117 
118  TPixel
120  {
121  return m_OutputMaximum;
122  }
123 
124  TPixel
125  operator()(TPixel input)
126  {
127  // Get the minimum and maximum output values
128  static const auto min = static_cast<float>(m_OutputMinimum);
129  static const auto max = static_cast<float>(m_OutputMaximum);
130 
131  // Compute the output
132  float output = static_cast<float>(input) + m_Mean + m_StandardDeviation * m_Generator->GetVariate();
133 
134  // Clamp the output value in valid range
135  output = (output < min ? min : output);
136  output = (output > max ? max : output);
137 
138  return static_cast<TPixel>(output);
139  }
140 
141 private:
144  float m_Mean;
147 };
148 
149 
170 template <class TInputImage>
171 class ITK_TEMPLATE_EXPORT AdditiveGaussianNoiseImageFilter : public itk::ImageToImageFilter<TInputImage, TInputImage>
172 {
173 public:
174  ITK_DISALLOW_COPY_AND_MOVE(AdditiveGaussianNoiseImageFilter);
175 
181 
183  itkNewMacro(Self);
184 
186 #ifdef itkOverrideGetNameOfClassMacro
187  itkOverrideGetNameOfClassMacro(AdditiveGaussianNoiseImageFilter);
188 #else
190 #endif
191 
192 
194  using OutputImageRegionType = typename Superclass::OutputImageRegionType;
195  using OutputImagePointer = typename Superclass::OutputImagePointer;
196 
198  using InputImageType = TInputImage;
199  using InputImagePointer = typename InputImageType::Pointer;
200  using InputImageConstPointer = typename InputImageType::ConstPointer;
201  using InputImageRegionType = typename InputImageType::RegionType;
202  using InputImagePixelType = typename InputImageType::PixelType;
203  using InputPixelType = typename InputImageType::PixelType;
204 
206  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
207 
208  // virtual void GenerateOutputInformation();
209 
210  void
211  GenerateData() override;
212 
213  // Accessor & Mutator methods
214 
219  void
220  SetMean(float mean)
221  {
222  m_NoiseFilter->GetFunctor().SetMean(mean);
223  this->Modified();
224  }
226 
231  float
232  GetMean() const
233  {
234  return m_NoiseFilter->GetFunctor().GetMean();
235  }
236 
241  void
242  SetStandardDeviation(float stddev)
243  {
244  m_NoiseFilter->GetFunctor().SetStandardDeviation(stddev);
245  this->Modified();
246  }
248 
253  float
255  {
256  return m_NoiseFilter->GetFunctor().GetStandardDeviation();
257  }
258 
265  void
266  SetSeed(unsigned long seed)
267  {
268  m_NoiseFilter->GetFunctor().SetSeed(seed);
269  this->Modified();
270  }
272 
274  void
276  {
277  if (min == m_NoiseFilter->GetFunctor().GetOutputMinimum())
278  {
279  return;
280  }
281  m_NoiseFilter->GetFunctor().SetOutputMinimum(min);
282  this->Modified();
283  }
285 
287  InputImagePixelType
289  {
290  return m_NoiseFilter->GetFunctor().GetOutputMinimum();
291  }
292 
294  void
296  {
297  if (max == m_NoiseFilter->GetFunctor().GetOutputMaximum())
298  {
299  return;
300  }
301  m_NoiseFilter->GetFunctor().SetOutputMaximum(max);
302  this->Modified();
303  }
305 
307  InputImagePixelType
309  {
310  return m_NoiseFilter->GetFunctor().GetOutputMaximum();
311  }
312 
313 protected:
315 
316  void
317  PrintSelf(std::ostream & os, itk::Indent indent) const override;
318 
319 public:
320  using NoiseFilterType = itk::UnaryFunctorImageFilter<InputImageType,
321  InputImageType,
323 
324 private:
326 };
327 
328 } /* end namespace rtk */
329 
330 #ifndef ITK_MANUAL_INSTANTIATION
331 # include "rtkAdditiveGaussianNoiseImageFilter.hxx"
332 #endif
333 
334 #endif /* rtkAdditiveGaussianNoiseImageFilter_h */
typename OutputImageType::Pointer OutputImagePointer
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::RegionType InputImageRegionType
Pixel functor that adds Gaussian noise.
typename OutputImageType::RegionType OutputImageRegionType
itk::Statistics::NormalVariateGenerator::Pointer m_Generator
typename InputImageType::PixelType InputImagePixelType