OTB  9.0.0
Orfeo Toolbox
otbAngularProjectionBinaryImageFilter.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 otbAngularProjectionBinaryImageFilter_hxx
22 #define otbAngularProjectionBinaryImageFilter_hxx
24 
25 #include <vnl/vnl_math.h>
26 
27 #include "itkImageRegionIterator.h"
28 #include "itkProgressReporter.h"
29 
30 namespace otb
31 {
32 
33 template <class TInputImage, class TOutputImage, class TPrecision>
35 {
36  this->SetNumberOfRequiredInputs(2);
37 }
38 
39 template <class TInputImage, class TOutputImage, class TPrecision>
41 {
42  this->SetNthInput(0, const_cast<InputImageType*>(inputPtr));
43 }
44 
45 template <class TInputImage, class TOutputImage, class TPrecision>
47 {
48  this->SetNthInput(1, const_cast<InputImageType*>(inputPtr));
49 }
50 
51 template <class TInputImage, class TOutputImage, class TPrecision>
53 {
54  if (this->GetNumberOfInputs() < 1)
55  {
56  return nullptr;
57  }
58 
59  return static_cast<const TInputImage*>(this->itk::ProcessObject::GetInput(0));
60 }
61 
62 template <class TInputImage, class TOutputImage, class TPrecision>
64 {
65  if (this->GetNumberOfInputs() < 2)
66  {
67  return nullptr;
68  }
69 
70  return static_cast<const TInputImage*>(this->itk::ProcessObject::GetInput(1));
71 }
72 
73 template <class TInputImage, class TOutputImage, class TPrecision>
75 {
76  m_AngleSet = angle;
77  this->SetNumberOfRequiredOutputs(angle.size());
78  for (unsigned int i = 0; i < this->GetNumberOfRequiredOutputs(); ++i)
79  {
80  this->SetNthOutput(i, OutputImageType::New());
81  }
82  this->Modified();
83 }
84 
85 template <class TInputImage, class TOutputImage, class TPrecision>
87 {
88  Superclass::GenerateOutputInformation();
89  for (unsigned int i = 0; i < this->GetNumberOfOutputs(); ++i)
90  {
91  this->GetOutput(i)->SetRegions(this->GetInput()->GetRequestedRegion());
92  }
93 }
94 
95 template <class TInputImage, class TOutputImage, class TPrecision>
97  itk::ThreadIdType threadId)
98 {
99  itk::ProgressReporter reporter(this, threadId, outputRegionForThread.GetNumberOfPixels());
100 
101  InputImageRegionType inputRegionForThread;
102  this->CallCopyOutputRegionToInputRegion(inputRegionForThread, outputRegionForThread);
103 
104  itk::ImageRegionConstIterator<InputImageType> iter1(this->GetInput1(), inputRegionForThread);
105  iter1.GoToBegin();
106 
107  itk::ImageRegionConstIterator<InputImageType> iter2(this->GetInput2(), inputRegionForThread);
108  iter2.GoToBegin();
109 
110  std::vector<itk::ImageRegionIterator<OutputImageType>> outIter(this->GetNumberOfOutputs());
111  for (unsigned int i = 0; i < outIter.size(); ++i)
112  {
113  outIter[i] = itk::ImageRegionIterator<OutputImageType>(this->GetOutput(i), outputRegionForThread);
114  outIter[i].GoToBegin();
115  }
116 
117  while (!iter1.IsAtEnd() && !iter2.IsAtEnd())
118  {
119  for (unsigned int i = 0; i < outIter.size(); ++i)
120  {
121  outIter[i].Set(std::cos(m_AngleSet[i]) * iter1.Get() + std::sin(m_AngleSet[i]) * iter2.Get());
122  ++outIter[i];
123  }
124 
125  ++iter1;
126  ++iter2;
127 
128  reporter.CompletedPixel();
129  }
130 }
131 
132 } // end of namespace otb
133 
134 #endif
otb::AngularProjectionBinaryImageFilter::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadID) override
Definition: otbAngularProjectionBinaryImageFilter.hxx:96
otb::AngularProjectionBinaryImageFilter::InputImageRegionType
InputImageType::RegionType InputImageRegionType
Definition: otbAngularProjectionBinaryImageFilter.h:59
otb::AngularProjectionBinaryImageFilter::SetInput1
void SetInput1(const InputImageType *)
Definition: otbAngularProjectionBinaryImageFilter.hxx:40
otb::AngularProjectionBinaryImageFilter::InputImageType
TInputImage InputImageType
Definition: otbAngularProjectionBinaryImageFilter.h:54
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::AngularProjectionBinaryImageFilter::GetInput2
const InputImageType * GetInput2() const
Definition: otbAngularProjectionBinaryImageFilter.hxx:63
otb::AngularProjectionBinaryImageFilter::SetInput2
void SetInput2(const InputImageType *)
Definition: otbAngularProjectionBinaryImageFilter.hxx:46
otb::AngularProjectionBinaryImageFilter::GenerateOutputInformation
void GenerateOutputInformation() override
Definition: otbAngularProjectionBinaryImageFilter.hxx:86
otb::AngularProjectionBinaryImageFilter::AngularProjectionBinaryImageFilter
AngularProjectionBinaryImageFilter()
Definition: otbAngularProjectionBinaryImageFilter.hxx:34
otb::AngularProjectionBinaryImageFilter::OutputImageRegionType
OutputImageType::RegionType OutputImageRegionType
Definition: otbAngularProjectionBinaryImageFilter.h:66
otb::AngularProjectionBinaryImageFilter::SetAngleSet
void SetAngleSet(std::vector< PrecisionType > &angle)
Definition: otbAngularProjectionBinaryImageFilter.hxx:74
otbAngularProjectionBinaryImageFilter.h
otb::AngularProjectionBinaryImageFilter::GetInput1
const InputImageType * GetInput1() const
Definition: otbAngularProjectionBinaryImageFilter.hxx:52