22 #ifndef __StreamingStatisticsMosaicFilter_H
23 #define __StreamingStatisticsMosaicFilter_H
29 #include "itkSimpleDataObjectDecorator.h"
30 #include "itkImageRegionConstIteratorWithOnlyIndex.h"
53 template <
class TInputImage,
class TOutputImage = TInputImage,
class TInternalValueType =
double>
85 typedef itk::ImageRegionConstIteratorWithOnlyIndex<OutputImageType>
IteratorType;
99 void AllocateOutputs()
override;
100 void ThreadedGenerateData(
const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
override;
101 void Reset()
override;
102 void Synthetize()
override;
107 using Superclass::MakeOutput;
112 return this->GetMeansOutput()->Get();
114 RealMatrixListObjectType* GetMeansOutput();
115 const RealMatrixListObjectType* GetMeansOutput()
const;
121 return this->GetStdsOutput()->Get();
123 RealMatrixListObjectType* GetStdsOutput();
124 const RealMatrixListObjectType* GetStdsOutput()
const;
130 return this->GetMinsOutput()->Get();
132 RealMatrixListObjectType* GetMinsOutput();
133 const RealMatrixListObjectType* GetMinsOutput()
const;
139 return this->GetMaxsOutput()->Get();
141 RealMatrixListObjectType* GetMaxsOutput();
142 const RealMatrixListObjectType* GetMaxsOutput()
const;
148 return this->GetMeansOfProductsOutput()->Get();
150 RealMatrixListObjectType* GetMeansOfProductsOutput();
151 const RealMatrixListObjectType* GetMeansOfProductsOutput()
const;
157 return this->GetAreasOutput()->Get();
159 RealMatrixObjectType* GetAreasOutput();
160 const RealMatrixObjectType* GetAreasOutput()
const;
188 Clear(nbOfBands, nbOfSamples);
203 void Clear(
unsigned int nbOfBands,
unsigned int nbOfSamples)
207 const InternalValueType infValue = itk::NumericTraits<InternalValueType>::NonpositiveMin();
220 unsigned int nbOfBands = pixel.Size();
223 for (
unsigned int band = 0; band < nbOfBands; band++)
229 if (pixelValue < m_min[band][sampleId])
230 m_min[band][sampleId] = pixelValue;
231 if (pixelValue > m_max[band][sampleId])
232 m_max[band][sampleId] = pixelValue;
235 m_sum[band][sampleId] += pixelValue;
236 m_sqSum[band][sampleId] += pixelValue * pixelValue;
243 Update(pixel_i, sampleId);
244 unsigned int nbOfBands = pixel_i.Size();
245 for (
unsigned int band = 0; band < nbOfBands; band++)
251 m_cosum[band][sampleId] += pixelValue_i * pixelValue_j;
258 unsigned int nbOfBands = other.
m_sum.rows();
259 unsigned int nbOfSamples = other.
m_sum.cols();
261 for (
unsigned int sampleId = 0; sampleId < nbOfSamples; sampleId++)
263 m_count[sampleId] += other.
m_count[sampleId];
264 for (
unsigned int band = 0; band < nbOfBands; band++)
266 m_sum[band][sampleId] += other.
m_sum[band][sampleId];
267 m_cosum[band][sampleId] += other.
m_cosum[band][sampleId];
268 m_sqSum[band][sampleId] += other.
m_sqSum[band][sampleId];
269 if (other.
m_min[band][sampleId] < m_min[band][sampleId])
270 m_min[band][sampleId] = other.
m_min[band][sampleId];
271 if (other.
m_max[band][sampleId] > m_max[band][sampleId])
272 m_max[band][sampleId] = other.
m_max[band][sampleId];
323 template <
class TInputImage,
class TOutputImage,
class TInternalValueType>
349 using Superclass::PushBackInput;
352 this->GetFilter()->PushBackInput(input);
358 return this->GetFilter()->GetMeansOutput()->Get();
362 return this->GetFilter()->GetMeansOutput();
366 return this->GetFilter()->GetMeansOutput();
373 return this->GetFilter()->GetStdsOutput()->Get();
377 return this->GetFilter()->GetStdsOutput();
381 return this->GetFilter()->GetStdsOutput();
388 return this->GetFilter()->GetMeansOfProductsOutput()->Get();
392 return this->GetFilter()->GetMeansOfProductsOutput();
396 return this->GetFilter()->GetMeansOfProductsOutput();
403 return this->GetFilter()->GetMinsOutput()->Get();
407 return this->GetFilter()->GetMinsOutput();
411 return this->GetFilter()->GetMinsOutput();
418 return this->GetFilter()->GetMaxsOutput()->Get();
422 return this->GetFilter()->GetMaxsOutput();
426 return this->GetFilter()->GetMaxsOutput();
433 return this->GetFilter()->GetAreasOutput()->Get();
437 return this->GetFilter()->GetAreasOutput();
441 return this->GetFilter()->GetAreasOutput();
461 #ifndef OTB_MANUAL_INSTANTIATION
This filter link a persistent filter with a StreamingImageVirtualWriter.
This filter is the base class for all mosaic filter persisting data through multiple update....
Superclass::OutputImageType OutputImageType
Superclass::InputImagePointer InputImagePointer
Superclass::InputImagePixelType InputImagePixelType
Superclass::InterpolatorPointerType InterpolatorPointerType
Superclass::OutputImageRegionType OutputImageRegionType
Superclass::InputImageInternalPixelType InputImageInternalPixelType
Superclass::InputImageType InputImageType
Superclass::InputImagePointType InputImagePointType
Superclass::OutputImagePointType OutputImagePointType
Superclass::InternalValueType InternalValueType
void Update(const ThreadResultsContainer &other)
void Clear(unsigned int nbOfBands, unsigned int nbOfSamples)
ThreadResultsContainer(const ThreadResultsContainer &other)
void Update(const InputImagePixelType &pixel_i, const InputImagePixelType &pixel_j, unsigned int sampleId)
void Update(const InputImagePixelType &pixel, unsigned int sampleId)
ThreadResultsContainer(unsigned int nbOfBands, unsigned int nbOfSamples)
Computes the statistics of a mosaic from an input images set. The output pixel value is equal to the ...
Superclass::InputImagePointer InputImagePointer
itk::ImageRegionConstIteratorWithOnlyIndex< OutputImageType > IteratorType
RealMatrixListType m_Maxs
itk::SmartPointer< Self > Pointer
RealMatrixListType GetStds() const
RealMatrixListType m_Mins
vnl_matrix< InternalValueType > RealMatrixType
itk::SmartPointer< const Self > ConstPointer
PersistentStatisticsMosaicFilter Self
itk::SimpleDataObjectDecorator< RealMatrixListType > RealMatrixListObjectType
virtual ~PersistentStatisticsMosaicFilter()
Superclass::InputImagePointType InputImagePointType
RealMatrixListType GetMaxs() const
Superclass::InputImageType InputImageType
RealMatrixListType GetMeansOfProducts() const
RealMatrixListType GetMins() const
Superclass::InterpolatorPointerType InterpolatorPointerType
Superclass::InternalValueType InternalValueType
Superclass::InputImageInternalPixelType InputImageInternalPixelType
RealMatrixType GetAreas() const
RealMatrixListType m_Means
PersistentStatisticsMosaicFilter(const Self &)
RealMatrixListType m_ProdMeans
Superclass::OutputImageRegionType OutputImageRegionType
Superclass::OutputImagePointType OutputImagePointType
RealMatrixListType GetMeans() const
PersistentMosaicFilter< TInputImage, TOutputImage, TInternalValueType > Superclass
itk::SimpleDataObjectDecorator< RealMatrixType > RealMatrixObjectType
std::vector< RealMatrixType > RealMatrixListType
RealMatrixListType m_Stds
void operator=(const Self &)
vnl_vector< InternalValueType > RealVectorType
itk::DataObject::Pointer DataObjectPointer
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
Superclass::InputImagePixelType InputImagePixelType
std::vector< ThreadResultsContainer > m_InternalThreadResults
Superclass::OutputImageType OutputImageType
This class streams the whole input image through the PersistentStatisticsMosaicFilter.
itk::SmartPointer< Self > Pointer
Superclass::FilterType::InputImageType InputImageType
RealMatrixListType GetMeans() const
Superclass::FilterType::RealMatrixType RealMatrixType
StreamingStatisticsMosaicFilter(const Self &)
RealMatrixListObjectType * GetMeansOutput()
const RealMatrixListObjectType * GetMaxsOutput() const
RealMatrixListObjectType * GetMaxsOutput()
~StreamingStatisticsMosaicFilter() override
RealMatrixListObjectType * GetStdsOutput()
RealMatrixType GetAreas() const
const RealMatrixListObjectType * GetMeansOfProductsOutput() const
RealMatrixListType GetStds() const
PersistentFilterStreamingDecorator< PersistentStatisticsMosaicFilter< TInputImage, TOutputImage, TInternalValueType > > Superclass
RealMatrixListObjectType * GetMeansOfProductsOutput()
StreamingStatisticsMosaicFilter Self
void operator=(const Self &)
const RealMatrixListObjectType * GetStdsOutput() const
RealMatrixListType GetMins() const
Superclass::FilterType::RealMatrixObjectType RealMatrixObjectType
const RealMatrixObjectType * GetAreasOutput() const
void PushBackInput(InputImageType *input)
const RealMatrixListObjectType * GetMeansOutput() const
itk::SmartPointer< const Self > ConstPointer
const RealMatrixListObjectType * GetMinsOutput() const
StreamingStatisticsMosaicFilter()
RealMatrixListType GetMaxs() const
Superclass::FilterType::RealMatrixListObjectType RealMatrixListObjectType
RealMatrixListType GetMeansOfProducts() const
RealMatrixObjectType * GetAreasOutput()
Superclass::FilterType::RealMatrixListType RealMatrixListType
RealMatrixListObjectType * GetMinsOutput()
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.