OTB  9.0.0
Orfeo Toolbox
otbStreamingShrinkImageFilter.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 otbStreamingShrinkImageFilter_hxx
22 #define otbStreamingShrinkImageFilter_hxx
23 
25 #include "itkImageRegionIteratorWithIndex.h"
26 #include "otbMacro.h"
27 #include "itkProgressReporter.h"
28 
29 namespace otb
30 {
31 
32 template <class TImage>
34 {
35 }
36 
37 template <class TImage>
39 {
40 }
41 
42 template <class TImage>
43 void StreamingShrinkStreamingManager<TImage>::PrepareStreaming(itk::DataObject* input, const RegionType& region)
44 {
45  typedef otb::StreamingShrinkImageRegionSplitter TileSplitterType;
46  TileSplitterType::Pointer splitter = TileSplitterType::New();
47  splitter->SetShrinkFactor(m_ShrinkFactor);
48  this->m_Splitter = splitter;
49 
50  unsigned long nbDivisions = this->EstimateOptimalNumberOfDivisions(input, region, 0);
51  this->m_ComputedNumberOfSplits = this->m_Splitter->GetNumberOfSplits(region, nbDivisions);
52  otbMsgDevMacro(<< "Number of split : " << this->m_ComputedNumberOfSplits)
53 
54  // Save the region to generate the splits later
55  this->m_Region = region;
56 }
57 
58 
60 template <class TInputImage, class TOutputImage>
62 {
63  this->SetNumberOfRequiredInputs(1);
64  this->SetNumberOfRequiredOutputs(1);
65 }
67 
69 template <class TInputImage, class TOutputImage>
71 {
72 }
73 
74 template <class TInputImage, class TOutputImage>
76 {
77  Superclass::GenerateOutputInformation();
78 
79  const InputImageType* input = this->GetInput();
80 
81  OutputImageType* output = this->GetOutput();
82 
83  if (input)
84  {
85  output->CopyInformation(input);
86  output->SetLargestPossibleRegion(input->GetLargestPossibleRegion());
87 
88  if (output->GetRequestedRegion().GetNumberOfPixels() == 0)
89  {
90  output->SetRequestedRegion(output->GetLargestPossibleRegion());
91  }
92  }
93 }
94 
95 template <class TInputImage, class TOutputImage>
97 {
98  // This is commented to prevent the streaming of the whole image for the first stream strip
99  // It shall not cause any problem because the output image of this filter is not intended to be used.
100  // InputImagePointer image = const_cast< TInputImage * >( this->GetInput() );
101  // this->GraftOutput( image );
102  // Nothing that needs to be allocated for the remaining outputs
103 }
104 
105 
106 template <class TInputImage, class TOutputImage>
108 {
109  // Get pointers to the input and output
110  InputImageType* inputPtr = const_cast<InputImageType*>(this->GetInput());
111  inputPtr->UpdateOutputInformation();
112 
113  m_ShrunkOutput = OutputImageType::New();
114  m_ShrunkOutput->CopyInformation(inputPtr);
115 
116  const typename InputImageType::SpacingType& inputSpacing = inputPtr->GetSignedSpacing();
117  const typename InputImageType::SizeType& inputSize = inputPtr->GetLargestPossibleRegion().GetSize();
118  const typename InputImageType::IndexType& inputIndex = inputPtr->GetLargestPossibleRegion().GetIndex();
119  typename InputImageType::IndexType startIndex;
120 
121  typename OutputImageType::SpacingType shrunkOutputSpacing;
122  typename OutputImageType::RegionType shrunkOutputLargestPossibleRegion;
123  typename OutputImageType::SizeType shrunkOutputSize;
124  typename OutputImageType::IndexType shrunkOutputStartIndex;
125  typename OutputImageType::PointType shrunkOutputOrigin;
126 
127  for (unsigned int i = 0; i < OutputImageType::ImageDimension; ++i)
128  {
129  startIndex[i] = inputIndex[i] + (m_ShrinkFactor - 1) / 2;
130  if (m_ShrinkFactor > inputSize[i])
131  startIndex[i] = inputIndex[i] + (inputSize[i] - 1) / 2;
132  m_Offset[i] = startIndex[i] % m_ShrinkFactor;
133  shrunkOutputSpacing[i] = inputSpacing[i] * static_cast<double>(m_ShrinkFactor);
134  shrunkOutputSize[i] = inputSize[i] > m_ShrinkFactor ? inputSize[i] / m_ShrinkFactor : 1;
135 
136  shrunkOutputOrigin[i] = inputPtr->GetOrigin()[i] + inputSpacing[i] * startIndex[i];
137 
138  // we choose to output a region with a start index [0,0]
139  // the origin is set accordingly
140  shrunkOutputStartIndex[i] = 0;
141  }
142 
143  m_ShrunkOutput->SetSignedSpacing(shrunkOutputSpacing);
144  m_ShrunkOutput->SetOrigin(shrunkOutputOrigin);
145 
146  shrunkOutputLargestPossibleRegion.SetSize(shrunkOutputSize);
147  shrunkOutputLargestPossibleRegion.SetIndex(shrunkOutputStartIndex);
148 
149  m_ShrunkOutput->SetRegions(shrunkOutputLargestPossibleRegion);
150  m_ShrunkOutput->Allocate();
151 }
152 
153 template <class TInputImage, class TOutputImage>
155 {
156 }
157 
158 template <class TInputImage, class TOutputImage>
160 {
161 }
162 
163 template <class TInputImage, class TOutputImage>
164 void PersistentShrinkImageFilter<TInputImage, TOutputImage>::ThreadedGenerateData(const RegionType& outputRegionForThread, itk::ThreadIdType threadId)
165 {
166  // std::cout << "outputRegionForThread " << threadId << " " << outputRegionForThread << std::endl;
167  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
168  const InputImageType* inputPtr = this->GetInput();
169 
170  itk::ImageRegionConstIteratorWithIndex<InputImageType> inIt(inputPtr, outputRegionForThread);
171  for (inIt.GoToBegin(); !inIt.IsAtEnd(); ++inIt, progress.CompletedPixel())
172  {
173  const IndexType& inIndex = inIt.GetIndex();
174  // TODO the pixel value should be taken near the centre of the cell, not at the corners
175  if ((inIndex[0] - m_Offset[0]) % m_ShrinkFactor == 0 && (inIndex[1] - m_Offset[1]) % m_ShrinkFactor == 0)
176  {
177  IndexType shrunkIndex;
178  shrunkIndex[0] = (inIndex[0] - m_Offset[0]) / m_ShrinkFactor;
179  shrunkIndex[1] = (inIndex[1] - m_Offset[1]) / m_ShrinkFactor;
180  if (m_ShrunkOutput->GetLargestPossibleRegion().IsInside(shrunkIndex))
181  m_ShrunkOutput->SetPixel(shrunkIndex, inIt.Get());
182  }
183  }
184 }
185 
186 template <class TInputImage, class TOutputImage>
188 {
189 }
190 
191 template <class TImage, class TOutputImage>
192 void PersistentShrinkImageFilter<TImage, TOutputImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
193 {
194  Superclass::PrintSelf(os, indent);
195  os << indent << "Shrink factor: " << m_ShrinkFactor << std::endl;
196 }
197 
198 } // End namespace otb
199 #endif
otb::PersistentShrinkImageFilter::Synthetize
void Synthetize(void) override
Definition: otbStreamingShrinkImageFilter.hxx:154
otb::PersistentShrinkImageFilter::Reset
void Reset(void) override
Definition: otbStreamingShrinkImageFilter.hxx:107
otb::StreamingShrinkImageRegionSplitter
Definition: otbStreamingShrinkImageFilter.h:34
otb::PersistentShrinkImageFilter::PersistentShrinkImageFilter
PersistentShrinkImageFilter()
Definition: otbStreamingShrinkImageFilter.hxx:61
otb::PersistentShrinkImageFilter::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbStreamingShrinkImageFilter.hxx:192
otb::PersistentShrinkImageFilter::ThreadedGenerateData
void ThreadedGenerateData(const RegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbStreamingShrinkImageFilter.hxx:164
otb::PersistentShrinkImageFilter::~PersistentShrinkImageFilter
~PersistentShrinkImageFilter() override
Definition: otbStreamingShrinkImageFilter.hxx:70
otb::PersistentShrinkImageFilter::AllocateOutputs
void AllocateOutputs() override
Definition: otbStreamingShrinkImageFilter.hxx:96
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbMacro.h
otb::StreamingShrinkStreamingManager::StreamingShrinkStreamingManager
StreamingShrinkStreamingManager()
Definition: otbStreamingShrinkImageFilter.hxx:33
otb::PersistentShrinkImageFilter::AfterThreadedGenerateData
void AfterThreadedGenerateData() override
Definition: otbStreamingShrinkImageFilter.hxx:187
otb::PersistentShrinkImageFilter::InputImageType
TInputImage InputImageType
Definition: otbStreamingShrinkImageFilter.h:183
otb::PersistentShrinkImageFilter::IndexType
TInputImage::IndexType IndexType
Definition: otbStreamingShrinkImageFilter.h:190
otb::StreamingShrinkStreamingManager::RegionType
ImageType::RegionType RegionType
Definition: otbStreamingShrinkImageFilter.h:126
otb::PersistentShrinkImageFilter::RegionType
TInputImage::RegionType RegionType
Definition: otbStreamingShrinkImageFilter.h:188
otb::PersistentShrinkImageFilter::GenerateOutputInformation
void GenerateOutputInformation() override
Definition: otbStreamingShrinkImageFilter.hxx:75
otb::StreamingShrinkStreamingManager::PrepareStreaming
void PrepareStreaming(itk::DataObject *input, const RegionType &region) override
Definition: otbStreamingShrinkImageFilter.hxx:43
otb::PersistentShrinkImageFilter::BeforeThreadedGenerateData
void BeforeThreadedGenerateData() override
Definition: otbStreamingShrinkImageFilter.hxx:159
otbMsgDevMacro
#define otbMsgDevMacro(x)
Definition: otbMacro.h:64
otb::StreamingShrinkStreamingManager::~StreamingShrinkStreamingManager
~StreamingShrinkStreamingManager() override
Definition: otbStreamingShrinkImageFilter.hxx:38
otb::PersistentShrinkImageFilter::OutputImageType
TOutputImage OutputImageType
Definition: otbStreamingShrinkImageFilter.h:194
otbStreamingShrinkImageFilter.h