OTB  9.0.0
Orfeo Toolbox
otbFastNLMeansImageFilter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2022 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbFastNLMeansImageFilter_hxx
22 #define otbFastNLMeansImageFilter_hxx
23 
25 #include "itkImageRegionIterator.h"
26 #include "itkImageRegionConstIterator.h"
27 #include "itkNumericTraits.h"
28 #include <vector>
29 #include <tuple>
30 
31 namespace otb
32 {
33  template<class TInputImage, class TOutputImage>
35  {
36  // Define default attributes values
37  m_HalfSearchSize.Fill(7);
38  m_HalfPatchSize.Fill(2);
39  m_Var = 0.;
40  m_CutoffDistance = 0.1;
41 
42  m_NormalizeDistance = m_CutoffDistance * m_CutoffDistance
43  * (2*m_HalfPatchSize[m_ROW]+1) * (2*m_HalfPatchSize[m_COL]+1);
44  }
45 
46  template<class TInputImage, class TOutputImage>
47  std::tuple< typename NLMeansFilter<TInputImage, TOutputImage>::InRegionType,
48  int, int, int, int, bool>
50  (const OutRegionType& outputRegion) const
51  {
52  InImageConstPointerType inputPtr = this->GetInput();
53  auto const& inputSize = inputPtr->GetLargestPossibleRegion().GetSize();
54 
55  // Get output region specification
56  auto const& outIndex = outputRegion.GetIndex();
57  auto const& outSize = outputRegion.GetSize();
58 
59  // Define margin for processing
60  const InSizeType halfMargin = m_HalfSearchSize + m_HalfPatchSize;
61  const InSizeType sizeTwo = {{2,2}};
62  const InSizeType fullMargin = sizeTwo*halfMargin;
63 
64  // Define region to read
65  InIndexType inIndex = outIndex - halfMargin;
66  InSizeType requestedSize = outSize + fullMargin;
67 
68  // Initialize parameters for mirror padding
69  bool needMirrorPadding = false;
70  int mirrorFirstRow = 0;
71  int mirrorFirstCol = 0;
72  int mirrorLastRow = 0;
73  int mirrorLastCol = 0;
74  // Check that the requested region is inside image boundaries
75  // If not, store number of missing data and update region
76  if (inIndex[m_COL] < 0){
77  needMirrorPadding = true;
78  mirrorFirstCol = -inIndex[m_COL];
79  inIndex[m_COL] = 0;
80  requestedSize[m_COL] -= mirrorFirstCol;
81  }
82  if (inIndex[m_ROW] < 0){
83  needMirrorPadding = true;
84  mirrorFirstRow = -inIndex[m_ROW];
85  inIndex[m_ROW] = 0;
86  requestedSize[m_ROW] -= mirrorFirstRow;
87  }
88  unsigned int lastCol = inIndex[m_COL] + requestedSize[m_COL];
89  if (lastCol >= inputSize[m_COL]){
90  needMirrorPadding = true;
91  mirrorLastCol = lastCol - inputSize[m_COL];
92  requestedSize[m_COL] -= mirrorLastCol;
93  }
94  unsigned int lastRow = inIndex[m_ROW] + requestedSize[m_ROW];
95  if (lastRow >= inputSize[m_ROW]){
96  needMirrorPadding = true;
97  mirrorLastRow = lastRow - inputSize[m_ROW];
98  requestedSize[m_ROW] -= mirrorLastRow;
99  }
100 
101  InRegionType inRequestedRegion(inIndex, requestedSize);
102  return std::make_tuple(inRequestedRegion, mirrorFirstRow, mirrorFirstCol,
103  mirrorLastRow, mirrorLastCol, needMirrorPadding);
104  }
105 
106  template <class TInputImage, class TOutputImage>
109  {
110  // Call the superclass' implementation of this method
111  Superclass::GenerateInputRequestedRegion();
112 
113  auto const& outputRequestedRegion = this->GetOutput()->GetRequestedRegion();
114  auto regionAndMirror = this->OutputRegionToInputRegion(outputRequestedRegion);
115  InRegionType inRequestedRegion = std::get<0>(regionAndMirror);
116 
117  InImageType * inputPtr = const_cast<InImageType * >(this->GetInput());
118  inputPtr->SetRequestedRegion(inRequestedRegion);
119  }
120 
121  template<class TInputImage, class TOutputImage>
122  void
124  (const OutRegionType& outputRegionForThread,
125  itk::ThreadIdType itkNotUsed(threadId))
126  {
127  InImageConstPointerType inputPtr = this->GetInput();
128  auto regionAndMirror = OutputRegionToInputRegion(outputRegionForThread);
129  // Unpack all values returned
130  InRegionType inputRegionForThread = std::get<0>(regionAndMirror);
131  int mirrorFirstRow = std::get<1>(regionAndMirror);
132  int mirrorFirstCol = std::get<2>(regionAndMirror);
133  int mirrorLastRow = std::get<3>(regionAndMirror);
134  int mirrorLastCol = std::get<4>(regionAndMirror);
135  bool needMirror = std::get<5>(regionAndMirror);
136 
137  // initialize and allocate vector to store temporary output values
138  // It makes it easier to store them in vectors to access various non-contiguous locations
139  auto const& outSize = outputRegionForThread.GetSize();
140  std::vector<double> outTemp(outSize[m_ROW]*outSize[m_COL]);
141  // initialize and allocate buffer to store all weights
142  std::vector<double> weights(outSize[m_ROW]*outSize[m_COL]);
143 
144  typedef itk::ImageRegionConstIterator<InImageType> InIteratorType;
145  InIteratorType inIt(inputPtr, inputRegionForThread);
146  auto const& inputSize = inputRegionForThread.GetSize();
147 
148  auto mirrorCol = inputSize[m_COL] + mirrorFirstCol + mirrorLastCol;
149  auto mirrorRow = inputSize[m_ROW] + mirrorFirstRow + mirrorLastRow;
150  InSizeType const& mirrorSize = {{mirrorCol, mirrorRow}};
151 
152  std::vector<double> dataInput(mirrorSize[m_ROW]*mirrorSize[m_COL]);
153  inIt.GoToBegin();
154  for (unsigned int row=static_cast<unsigned int>(mirrorFirstRow);
155  row<static_cast<unsigned int>(mirrorFirstRow)+inputSize[m_ROW]; row++)
156  for (unsigned int col=static_cast<unsigned int>(mirrorFirstCol);
157  col<static_cast<unsigned int>(mirrorFirstCol)+inputSize[m_COL]; col++)
158  {
159  auto index = row * mirrorSize[m_COL] + col;
160  dataInput[index] = static_cast<double>(inIt.Get());
161  ++inIt;
162  }
163 
164  if (needMirror)
165  {
166  // Perform mirror on upper lines
167  for (int row=0; row<mirrorFirstRow; row++)
168  {
169  int lineToCopy = (2*mirrorFirstRow - row)*mirrorSize[m_COL];
170  std::copy(dataInput.begin() + lineToCopy,
171  dataInput.begin() + lineToCopy + mirrorSize[m_COL],
172  dataInput.begin() + row*mirrorSize[m_COL] );
173  }
174  // Perform mirror on lower lines
175  int lastRowRead = mirrorFirstRow+inputSize[m_ROW];
176  for (int row=0; row<mirrorLastRow; row++)
177  {
178  int lineToCopy = (lastRowRead - row -2)*mirrorSize[m_COL];
179  std::copy(dataInput.begin() + lineToCopy,
180  dataInput.begin() + lineToCopy + mirrorSize[m_COL],
181  dataInput.begin() + (lastRowRead + row)*mirrorSize[m_COL]);
182  }
183  // Perform mirror on left-hand columns
184  if (mirrorFirstCol > 0) {
185  for (unsigned int row=0; row<mirrorSize[m_ROW]; row++)
186  {
187  std::reverse_copy(dataInput.begin() + row*mirrorSize[m_COL] + mirrorFirstCol+1,
188  dataInput.begin() + row*mirrorSize[m_COL] +2*mirrorFirstCol+1,
189  dataInput.begin() + row*mirrorSize[m_COL]);
190  }
191 
192  }
193  // Perform mirror on right-hand columns
194  if (mirrorLastCol > 0){
195  for (unsigned int row=0; row<mirrorSize[m_ROW]; row++)
196  {
197  std::reverse_copy(dataInput.begin() + (row+1)*mirrorSize[m_COL] - 2*mirrorLastCol-1,
198  dataInput.begin() + (row+1)*mirrorSize[m_COL] - mirrorLastCol-1,
199  dataInput.begin() + (row+1)*mirrorSize[m_COL] - mirrorLastCol);
200  }
201  }
202  }
203 
204  // For loops on all shifts possible
205  int fullMarginRow = static_cast<int>(m_HalfSearchSize[m_ROW]+m_HalfPatchSize[m_ROW]);
206  int fullMarginCol = static_cast<int>(m_HalfSearchSize[m_COL]+m_HalfPatchSize[m_COL]);
207  int searchSizeRow = static_cast<int>(m_HalfSearchSize[m_ROW]);
208  int searchSizeCol = static_cast<int>(m_HalfSearchSize[m_COL]);
209  // Allocate integral image
210  const InSizeType sizeTwo = {{2,2}};
211  auto const& inSize = outSize + sizeTwo * m_HalfPatchSize;
212 
213  std::vector<double> imIntegral(inSize[m_ROW]*inSize[m_COL]);
214  for (int drow=-searchSizeRow; drow < searchSizeRow+1; drow++)
215  for (int dcol=-searchSizeCol; dcol < searchSizeCol+1; dcol++)
216  {
217  // Compute integral image for current shift (drow, dcol)
218  OutIndexType shift = {{dcol, drow}};
219  ComputeIntegralImage(dataInput, imIntegral, shift, inSize, mirrorSize);
220 
221  for(unsigned int row=0; row<outSize[m_ROW]; row++)
222  for (unsigned int col=0; col<outSize[m_COL]; col++)
223  {
224  // Compute distance from integral image for patch centered at
225  // (row, col) + (m_HalfPatchSize, m_HalfPatchSize)
226  OutPixelType distance = ComputeDistance(row, col, imIntegral, inSize[m_COL]);
227  if (distance < 5.0)
228  {
229  double weight = exp(static_cast<double>(-distance));
230  outTemp[row*outSize[m_COL] + col] += weight*dataInput[(row+drow+fullMarginRow)*mirrorSize[m_COL]
231  + col+dcol+fullMarginCol];
232  weights[row*outSize[m_COL] + col] += weight;
233  }
234  }
235  }
236 
237  // Normalize all results by dividing output by weights (store in output)
238  typedef itk::ImageRegionIterator<OutImageType> OutputIteratorType;
239  OutImagePointerType outputPtr = this->GetOutput();
240  OutputIteratorType outIt(outputPtr, outputRegionForThread);
241  outIt.GoToBegin();
242  for(unsigned int index=0; index<outSize[m_ROW]*outSize[m_COL]; index++)
243  {
244  outIt.Set(static_cast<OutPixelType>(outTemp[index]/weights[index]));
245  ++outIt;
246  }
247  }
248 
249  template<class TInputImage, class TOutputImage>
250  void
252  (const std::vector<double> & dataInput,
253  std::vector<double> &imIntegral,
254  const OutIndexType shift, const InSizeType sizeIntegral, const InSizeType sizeInput) const
255  {
256  // dataInput has a margin of m_HalfSearchSize+m_HalfPatchSize to allow
257  // computation of all shifts (computation of all integral images)
258  // integral images just have the m_HalfPatchSize margin necessary
259  // to compute patches differences for a given shift
260  // hence, the first point used in computation for the non-shifted image
261  // is located at m_HalfSearchSize
262  auto const& offsetRef = m_HalfSearchSize;
263  OutSizeType const& offsetShift = {{offsetRef[0] + shift[0], offsetRef[1] + shift[1]}};
264 
265  // Initialize integral image (compute position (0,0))
266  auto indexInput = offsetRef[m_ROW]*sizeInput[m_COL] + offsetRef[m_COL];
267  auto indexShift = offsetShift[m_ROW]*sizeInput[m_COL] + offsetShift[m_COL];
268  double diff = dataInput[indexInput] - dataInput[indexShift];
269  imIntegral[0] = (diff * diff) - m_Var;
270  // Compute first line of integral image
271  for (unsigned int col=1; col<sizeIntegral[m_COL]; col++)
272  {
273  auto indexInputCol = indexInput + col;
274  auto indexShiftCol = indexShift + col;
275  diff = dataInput[indexInputCol] - dataInput[indexShiftCol];
276  double distance = diff * diff - m_Var;
277  imIntegral[col] = distance + imIntegral[col-1];
278  assert(imIntegral[col] < itk::NumericTraits<double>::max());
279  }
280  // Compute first column of integral image
281  for (unsigned int row=1; row<sizeIntegral[m_ROW]; row++)
282  {
283  auto indexInputRow = indexInput + row*sizeInput[m_COL];
284  auto indexShiftRow = indexShift + row*sizeInput[m_COL];
285  diff = dataInput[indexInputRow] - dataInput[indexShiftRow];
286  double distance = diff * diff - m_Var;
287  imIntegral[row*sizeIntegral[m_COL]] = distance + imIntegral[(row-1)*sizeIntegral[m_COL]];
288  assert(imIntegral[row*sizeIntegral[m_COL]] < itk::NumericTraits<double>::max());
289  }
290  // All initializations have been done previously
291  // Remaining coefficients can be computed
292  for (unsigned int row=1; row<sizeIntegral[m_ROW]; row++)
293  for (unsigned int col=1; col<sizeIntegral[m_COL]; col++)
294  {
295  indexInput = (offsetRef[m_ROW]+row)*sizeInput[m_COL] + offsetRef[m_COL]+col;
296  indexShift = (offsetShift[m_ROW]+row)*sizeInput[m_COL] + offsetShift[m_COL]+col;
297  diff = dataInput[indexInput] - dataInput[indexShift];
298  double distance = diff*diff - m_Var;
299  imIntegral[row*sizeIntegral[m_COL] + col] = distance + imIntegral[row*sizeIntegral[m_COL] + col-1]
300  + imIntegral[(row-1)*sizeIntegral[m_COL] + col] - imIntegral[(row-1)*sizeIntegral[m_COL] + col-1];
301  assert(imIntegral[row*sizeIntegral[m_COL] + col] < itk::NumericTraits<double>::max());
302  }
303 
304  }
305 
306  template <class TInputImage, class TOutputImage>
308  NLMeansFilter<TInputImage, TOutputImage>::ComputeDistance
309  (const unsigned int row, const unsigned int col,
310  const std::vector<double>& imIntegral, const unsigned int nbCols) const
311  {
312  // (row, col) is the central position of the local window in the output image
313  // however, integral image is shifted by (m_HalfPatchSize, m_HalfPatchSize) compared to output image
314  // Thus, (row, col) corresponds, in integral image, to the upper left corner of the local window
315  double distance_patch =
316  imIntegral[(row+2*m_HalfPatchSize[m_ROW])*nbCols + col+2*m_HalfPatchSize[m_COL]]
317  - imIntegral[row*nbCols + col+2*m_HalfPatchSize[m_COL]]
318  - imIntegral[(row+2*m_HalfPatchSize[m_ROW])*nbCols + col]
319  + imIntegral[row*nbCols + col];
320 
321  distance_patch = std::max(distance_patch, 0.0) / (m_NormalizeDistance);
322  return static_cast<OutPixelType>(distance_patch);
323  }
324 
325  template<class TInputImage, class TOutputImage>
327  ::PrintSelf(std::ostream & os, itk::Indent indent) const
328  {
329  Superclass::PrintSelf(os, indent);
330 
331  os<<indent<<"NL Means Patch Size : "<<2*m_HalfPatchSize[m_ROW]+1
332  <<" x "<<2*m_HalfPatchSize[m_COL]+1<< std::endl;
333  os<<indent<<"NL Means Window Search Size : "<<2*m_HalfSearchSize[m_ROW]+1
334  <<" x "<<2*m_HalfSearchSize[m_COL]+1<< std::endl;
335  os<<indent<<"NL Means variance : "<<m_Var<<std::endl;
336  os<<indent<<"NL Means threshold for similarity : "<<m_CutoffDistance
337  << std::endl;
338  }
339 
340 } // end namespace otb
341 
342 #endif
otb::NLMeansFilter::OutImagePointerType
OutImageType::Pointer OutImagePointerType
Definition: otbFastNLMeansImageFilter.h:68
otb::NLMeansFilter
This class implements a fast approximated version of NL Means denoising algorithm....
Definition: otbFastNLMeansImageFilter.h:48
otbFastNLMeansImageFilter.h
otb::NLMeansFilter::InIndexType
InImageType::IndexType InIndexType
Definition: otbFastNLMeansImageFilter.h:65
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::NLMeansFilter::OutRegionType
OutImageType::RegionType OutRegionType
Definition: otbFastNLMeansImageFilter.h:69
otb::NLMeansFilter::InRegionType
InImageType::RegionType InRegionType
Definition: otbFastNLMeansImageFilter.h:64
otb::NLMeansFilter::NLMeansFilter
NLMeansFilter()
Definition: otbFastNLMeansImageFilter.hxx:34
otb::NLMeansFilter::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbFastNLMeansImageFilter.hxx:327
otb::NLMeansFilter::InImageType
TInputImage InImageType
Definition: otbFastNLMeansImageFilter.h:58
otb::NLMeansFilter::GenerateInputRequestedRegion
void GenerateInputRequestedRegion() override
Definition: otbFastNLMeansImageFilter.hxx:108
otb::NLMeansFilter::InSizeType
InImageType::SizeType InSizeType
Definition: otbFastNLMeansImageFilter.h:66
otb::NLMeansFilter::OutIndexType
OutImageType::IndexType OutIndexType
Definition: otbFastNLMeansImageFilter.h:72
otb::NLMeansFilter::OutputRegionToInputRegion
std::tuple< InRegionType, int, int, int, int, bool > OutputRegionToInputRegion(const OutRegionType &outputRegion) const
Definition: otbFastNLMeansImageFilter.hxx:50
otb::NLMeansFilter::ThreadedGenerateData
void ThreadedGenerateData(const OutRegionType &outputRegionForThread, itk::ThreadIdType) override
Definition: otbFastNLMeansImageFilter.hxx:124
otb::NLMeansFilter::InImageConstPointerType
InImageType::ConstPointer InImageConstPointerType
Definition: otbFastNLMeansImageFilter.h:63
otb::NLMeansFilter::OutPixelType
OutImageType::PixelType OutPixelType
Definition: otbFastNLMeansImageFilter.h:70