1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289
|
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
Some parts of this code are derived from ITK. See ITKCopyright.txt
for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#ifndef otbStreamingMatrixTransposeMatrixImageFilter_h
#define otbStreamingMatrixTransposeMatrixImageFilter_h
#include "otbPersistentImageFilter.h"
#include "itkSimpleDataObjectDecorator.h"
#include "otbStreamingTraits.h"
#include "itkVariableSizeMatrix.h"
#include "itkVariableLengthVector.h"
#include "otbPersistentFilterStreamingDecorator.h"
namespace otb
{
/** \class PersistentMatrixTransposeMatrixImageFilter
* \brief Compute \f[X^T.Y \f]. Allow a padding of ones.
*
* \f[X\f] and \f[Y\f] are the input images.
* The padding has the effect of adding a component filled with ones to the image
*
* This filter persists its temporary data. It means that if you Update it n times on n different
* requested regions, the output statistics will be the result for the whole set of n regions.
*
* To reset the temporary data, one should call the Reset() function.
*
* To get the statistics once the regions have been processed via the pipeline, use the Synthetize() method.
*
* \sa StreamingTraits
* \sa StatisticsImageFilter
* \ingroup MathematicalStatisticsImageFilters
*
* \ingroup OTBImageManipulation
*/
template<class TInputImage, class TInputImage2>
class ITK_EXPORT PersistentMatrixTransposeMatrixImageFilter :
public PersistentImageFilter<TInputImage, TInputImage>
{
public:
/** Standard Self typedef */
typedef PersistentMatrixTransposeMatrixImageFilter Self;
typedef PersistentImageFilter<TInputImage, TInputImage> Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** Runtime information support. */
itkTypeMacro(PersistentMatrixTransposeMatrixImageFilter, PersistentImageFilter);
/** Image related typedefs. */
// First Input
typedef TInputImage ImageType;
typedef typename TInputImage::Pointer InputImagePointer;
typedef typename TInputImage::RegionType RegionType;
typedef typename TInputImage::SizeType SizeType;
typedef typename TInputImage::IndexType IndexType;
typedef typename TInputImage::PixelType PixelType;
typedef typename TInputImage::InternalPixelType InternalPixelType;
typedef typename TInputImage2::IndexType IndexType2;
typedef typename TInputImage2::PixelType PixelType2;
typedef typename TInputImage2::InternalPixelType InternalPixelType2;
itkStaticConstMacro(InputImageDimension, unsigned int, TInputImage::ImageDimension);
itkSetMacro(UsePadFirstInput, bool);
itkGetMacro(UsePadFirstInput, bool);
itkSetMacro(UsePadSecondInput, bool);
itkGetMacro(UsePadSecondInput, bool);
/** Image related typedefs. */
itkStaticConstMacro(ImageDimension, unsigned int, TInputImage::ImageDimension);
/** Type to use for computations. */
// First Input
typedef double RealType;
typedef itk::VariableLengthVector<RealType> RealPixelType;
/** Smart Pointer type to a DataObject. */
typedef typename itk::DataObject::Pointer DataObjectPointer;
typedef itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType;
/** Type of DataObjects used for scalar outputs */
typedef typename itk::Array<long> ArrayLongPixelType;
typedef typename itk::VariableSizeMatrix<RealType> MatrixType;
typedef typename std::vector<MatrixType> ArrayMatrixType;
typedef typename std::vector<RealPixelType> ArrayRealPixelType;
typedef typename std::vector<PixelType> ArrayPixelType;
typedef itk::SimpleDataObjectDecorator<RealPixelType> RealPixelObjectType;
typedef itk::SimpleDataObjectDecorator<PixelType> PixelObjectType;
typedef itk::SimpleDataObjectDecorator<MatrixType> MatrixObjectType;
/** Return the computed transpose(Image1)*Image2. */
MatrixType GetResult() const
{
return this->GetResultOutput()->Get();
}
MatrixObjectType* GetResultOutput();
const MatrixObjectType* GetResultOutput() const;
/** Make a DataObject of the correct type to be used as the specified
* output.
*/
DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) ITK_OVERRIDE;
using Superclass::MakeOutput;
/** Pass the input through unmodified. Do this by Grafting in the
* AllocateOutputs method.
*/
void AllocateOutputs() ITK_OVERRIDE;
void GenerateOutputInformation() ITK_OVERRIDE;
// Override since the filter needs all the data for the algorithm
void GenerateInputRequestedRegion() ITK_OVERRIDE;
void Reset(void) ITK_OVERRIDE;
void Synthetize(void) ITK_OVERRIDE;
/** Input wrapper */
void SetFirstInput(const TInputImage *input1)
{
this->SetInput(0, input1);
}
void SetSecondInput(const TInputImage2 *input2)
{
this->SetInput(1, input2);
}
const TInputImage* GetFirstInput()
{
if (this->GetNumberOfInputs() < 1)
{
return ITK_NULLPTR;
}
else return (static_cast<const TInputImage *>(this->itk::ProcessObject::GetInput(0)));
}
const TInputImage2* GetSecondInput()
{
if (this->GetNumberOfInputs() < 2)
{
return ITK_NULLPTR;
}
else return (static_cast<const TInputImage2 *>(this->itk::ProcessObject::GetInput(1)));
}
protected:
PersistentMatrixTransposeMatrixImageFilter();
~PersistentMatrixTransposeMatrixImageFilter() ITK_OVERRIDE {}
void PrintSelf(std::ostream& os, itk::Indent indent) const ITK_OVERRIDE;
/** Multi-thread version GenerateData. */
void ThreadedGenerateData(const RegionType& outputRegionForThread, itk::ThreadIdType threadId) ITK_OVERRIDE;
private:
PersistentMatrixTransposeMatrixImageFilter(const Self &); //purposely not implemented
void operator =(const Self&); //purposely not implemented
ArrayMatrixType m_ThreadSum;
bool m_UsePadFirstInput;
bool m_UsePadSecondInput;
/** Nulber Of Component per Pixel. Change for padding */
unsigned int m_NumberOfComponents1;
unsigned int m_NumberOfComponents2;
}; // end of class
/**===========================================================================*/
/** \class StreamingMatrixTransposeMatrixImageFilter
* \brief This class streams the whole input image through the PersistentMatrixTransposeMatrixImageFilter.
*
* This way, it allows computing \f[X^T.Y \f] where \f[X\f] and \f[Y\f] are the input images.
* first order global statistics of this image. It calls the Reset() method of the
* PersistentMatrixTransposeMatrixImageFilter before streaming the image and the
* Synthetize() method of the PersistentMatrixTransposeMatrixImageFilter after having streamed the image
* to compute the statistics. The accessor on the results are wrapping the accessors of the
* internal PersistentStatisticsImageFilter. The accessor on the pad options are also provided.
*
* \sa PersistentMatrixTransposeMatrixImageFilter
* \sa PersistentImageFilter
* \sa PersistentFilterStreamingDecorator
* \sa StreamingImageVirtualWriter
* \ingroup Streamed
* \ingroup Multithreaded
* \ingroup MathematicalStatisticsImageFilters
*
* \ingroup OTBImageManipulation
*/
template<class TInputImage1, class TInputImage2>
class ITK_EXPORT StreamingMatrixTransposeMatrixImageFilter :
public PersistentFilterStreamingDecorator
<PersistentMatrixTransposeMatrixImageFilter<TInputImage1, TInputImage2> >
{
public:
/** Standard Self typedef */
typedef StreamingMatrixTransposeMatrixImageFilter Self;
typedef PersistentFilterStreamingDecorator
<PersistentMatrixTransposeMatrixImageFilter<TInputImage1, TInputImage2> > Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Type macro */
itkNewMacro(Self);
/** Creation through object factory macro */
itkTypeMacro(StreamingMatrixTransposeMatrixImageFilter, PersistentFilterStreamingDecorator);
typedef typename Superclass::FilterType MatrixToTransposeMatrixFilterType;
typedef typename MatrixToTransposeMatrixFilterType::MatrixType MatrixType;
typedef itk::SimpleDataObjectDecorator<MatrixType> MatrixObjectType;
typedef TInputImage1 InputImage1Type;
typedef TInputImage2 InputImage2Type;
/** Interfaces to the embedded filter */
void SetFirstInput(InputImage1Type * input1)
{
this->GetFilter()->SetFirstInput(input1);
}
void SetSecondInput(InputImage2Type * input2)
{
this->GetFilter()->SetSecondInput(input2);
}
void SetUsePadFirstInput(bool pad)
{
this->GetFilter()->SetUsePadFirstInput(pad);
}
void SetUsePadSecondInput(bool pad)
{
this->GetFilter()->SetUsePadSecondInput(pad);
}
bool GetUsePadFirstInput(void)
{
return this->GetFilter()->GetUsePadFirstInput();
}
bool GetUsePadSecondInput(void)
{
return this->GetFilter()->GetUsePadSecondInput();
}
/** Return the computed transpose(Image1)*Image2. */
MatrixType GetResult(void) const
{
return this->GetResultOutput()->Get();
}
MatrixObjectType* GetResultOutput(void)
{
return this->GetFilter()->GetResultOutput();
}
const MatrixObjectType* GetResultOutput() const
{
return this->GetFilter()->GetResultOutput();
}
protected:
/** Constructor */
StreamingMatrixTransposeMatrixImageFilter() {};
/** Destructor */
~StreamingMatrixTransposeMatrixImageFilter() ITK_OVERRIDE {}
private:
StreamingMatrixTransposeMatrixImageFilter(const Self &); //purposely not implemented
void operator =(const Self&); //purposely not implemented
};
} // end namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbStreamingMatrixTransposeMatrixImageFilter.txx"
#endif
#endif
|