OTB  9.0.0
Orfeo Toolbox
List of all members
otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision > Class Template Reference

#include <otbConvolutionImageFilter.h>

+ Inheritance diagram for otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >:
+ Collaboration diagram for otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >:
typedef TInputImage InputImageType
 
typedef TOutputImage OutputImageType
 
typedef ConvolutionImageFilter Self
 
typedef itk::ImageToImageFilter< InputImageType, OutputImageTypeSuperclass
 
typedef itk::SmartPointer< SelfPointer
 
typedef itk::SmartPointer< const SelfConstPointer
 
typedef InputImageType::PixelType InputPixelType
 
typedef OutputImageType::PixelType OutputPixelType
 
typedef itk::NumericTraits< InputPixelType >::RealType InputRealType
 
typedef InputImageType::RegionType InputImageRegionType
 
typedef OutputImageType::RegionType OutputImageRegionType
 
typedef InputImageType::SizeType InputSizeType
 
typedef TFilterPrecision FilterPrecisionType
 
typedef itk::Array< FilterPrecisionTypeArrayType
 
typedef TBoundaryCondition BoundaryConditionType
 
static const unsigned int InputImageDimension = TInputImage::ImageDimension
 
static const unsigned int OutputImageDimension = TOutputImage::ImageDimension
 
InputSizeType m_Radius
 
ArrayType m_Filter
 
bool m_NormalizeFilter
 
static Pointer New ()
 
virtual ::itk::LightObject::Pointer CreateAnother (void) const
 
virtual const char * GetNameOfClass () const
 
virtual void SetRadius (const InputSizeType rad)
 
virtual const InputSizeTypeGetRadius () const
 
virtual void SetFilter (ArrayType filter)
 
virtual const ArrayTypeGetFilter () const
 
virtual void SetNormalizeFilter (bool _arg)
 
virtual bool GetNormalizeFilter ()
 
virtual void NormalizeFilterOn ()
 
virtual void NormalizeFilterOff ()
 
 typedef (itk::Concept::HasNumericTraits< InputPixelType >) InputHasNumericTraitsCheck
 
 ConvolutionImageFilter ()
 
 ~ConvolutionImageFilter () override
 
void PrintSelf (std::ostream &os, itk::Indent indent) const override
 
void ThreadedGenerateData (const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
 
void GenerateInputRequestedRegion () override
 
 ConvolutionImageFilter (const Self &)=delete
 
void operator= (const Self &)=delete
 

Detailed Description

template<class TInputImage, class TOutputImage, class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
class otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >

Applies a convolution filter to a mono channel image.

Computes an image which is the convolution of the input image with a filter.

The radius of the input filter is provided by the

SetInput()

method and the filters coefficients are given by an itk::Array passed to the

method.

By default, the input filter is not normalized but it can be using the NormalizeFilterOn() method.

This filter allows the user to choose the boundary conditions in the template parameters. Default boundary conditions are zero flux Neumann boundary conditions.

An optimized version of this filter using FFTW is available in the Orfeo ToolBox and will significantly improves performances especially for large kernels (see OverlapSaveConvolutionImageFilter).

See also
Image
Neighborhood
NeighborhoodOperator
NeighborhoodIterator
ImageBoundaryCondition
ZeroFluxNeumannBoundaryCondition
OverlapSaveConvolutionImageFilter

Definition at line 68 of file otbConvolutionImageFilter.h.

Member Typedef Documentation

◆ ArrayType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef itk::Array<FilterPrecisionType> otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::ArrayType

Convenient typedefs for simplifying declarations.

Definition at line 100 of file otbConvolutionImageFilter.h.

◆ BoundaryConditionType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef TBoundaryCondition otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::BoundaryConditionType

Convenient typedefs for simplifying declarations.

Definition at line 101 of file otbConvolutionImageFilter.h.

◆ ConstPointer

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef itk::SmartPointer<const Self> otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::ConstPointer

Convenient typedefs for simplifying declarations.

Definition at line 84 of file otbConvolutionImageFilter.h.

◆ FilterPrecisionType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef TFilterPrecision otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::FilterPrecisionType

Convenient typedefs for simplifying declarations.

Definition at line 99 of file otbConvolutionImageFilter.h.

◆ InputImageRegionType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef InputImageType::RegionType otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::InputImageRegionType

Convenient typedefs for simplifying declarations.

Definition at line 96 of file otbConvolutionImageFilter.h.

◆ InputImageType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef TInputImage otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::InputImageType

Convenient typedefs for simplifying declarations.

Definition at line 77 of file otbConvolutionImageFilter.h.

◆ InputPixelType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef InputImageType::PixelType otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::InputPixelType

Image typedef support.

Definition at line 90 of file otbConvolutionImageFilter.h.

◆ InputRealType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef itk::NumericTraits<InputPixelType>::RealType otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::InputRealType

Convenient typedefs for simplifying declarations.

Definition at line 95 of file otbConvolutionImageFilter.h.

◆ InputSizeType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef InputImageType::SizeType otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::InputSizeType

Convenient typedefs for simplifying declarations.

Definition at line 98 of file otbConvolutionImageFilter.h.

◆ OutputImageRegionType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef OutputImageType::RegionType otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::OutputImageRegionType

Convenient typedefs for simplifying declarations.

Definition at line 97 of file otbConvolutionImageFilter.h.

◆ OutputImageType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef TOutputImage otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::OutputImageType

Convenient typedefs for simplifying declarations.

Definition at line 78 of file otbConvolutionImageFilter.h.

◆ OutputPixelType

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef OutputImageType::PixelType otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::OutputPixelType

Convenient typedefs for simplifying declarations.

Definition at line 94 of file otbConvolutionImageFilter.h.

◆ Pointer

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef itk::SmartPointer<Self> otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::Pointer

Convenient typedefs for simplifying declarations.

Definition at line 83 of file otbConvolutionImageFilter.h.

◆ Self

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef ConvolutionImageFilter otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::Self

Standard class typedefs.

Definition at line 81 of file otbConvolutionImageFilter.h.

◆ Superclass

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
typedef itk::ImageToImageFilter<InputImageType, OutputImageType> otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::Superclass

Convenient typedefs for simplifying declarations.

Definition at line 82 of file otbConvolutionImageFilter.h.

Constructor & Destructor Documentation

◆ ConvolutionImageFilter() [1/2]

template<class TInputImage , class TOutputImage , class TBoundaryCondition , class TFilterPrecision >
otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::ConvolutionImageFilter
protected

End concept checking

Definition at line 38 of file otbConvolutionImageFilter.hxx.

◆ ~ConvolutionImageFilter()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::~ConvolutionImageFilter ( )
inlineoverrideprotected

Convenient typedefs for simplifying declarations.

Definition at line 158 of file otbConvolutionImageFilter.h.

◆ ConvolutionImageFilter() [2/2]

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::ConvolutionImageFilter ( const Self )
privatedelete

Convenient typedefs for simplifying declarations.

Member Function Documentation

◆ CreateAnother()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual::itk::LightObject::Pointer otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::CreateAnother ( void  ) const

Convenient typedefs for simplifying declarations.

◆ GenerateInputRequestedRegion()

template<class TInputImage , class TOutputImage , class TBoundaryCondition , class TFilterPrecision >
void otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::GenerateInputRequestedRegion
overrideprotected

ConvolutionImageFilter needs a larger input requested region than the output requested region. As such, ConvolutionImageFilter needs to provide an implementation for GenerateInputRequestedRegion() in order to inform the pipeline execution model.

See also
ImageToImageFilter::GenerateInputRequestedRegion()

Definition at line 47 of file otbConvolutionImageFilter.hxx.

References otbMsgDevMacro.

◆ GetFilter()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual const ArrayType& otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::GetFilter ( ) const
virtual

Convenient typedefs for simplifying declarations.

◆ GetNameOfClass()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual const char* otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

◆ GetNormalizeFilter()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual bool otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::GetNormalizeFilter ( )
virtual

Convenient typedefs for simplifying declarations.

◆ GetRadius()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual const InputSizeType& otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::GetRadius ( ) const
virtual

Get the radius of the neighborhood of the filter

◆ New()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
static Pointer otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::New ( )
static

Method for creation through the object factory.

◆ NormalizeFilterOff()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual void otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::NormalizeFilterOff ( )
virtual

Convenient typedefs for simplifying declarations.

◆ NormalizeFilterOn()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual void otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::NormalizeFilterOn ( )
virtual

Convenient typedefs for simplifying declarations.

◆ operator=()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
void otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::operator= ( const Self )
privatedelete

Convenient typedefs for simplifying declarations.

◆ PrintSelf()

template<class TInputImage , class TOutput , class TBoundaryCondition , class TFilterPrecision >
void otb::ConvolutionImageFilter< TInputImage, TOutput, TBoundaryCondition, TFilterPrecision >::PrintSelf ( std::ostream &  os,
itk::Indent  indent 
) const
overrideprotected

Standard "PrintSelf" method

Definition at line 168 of file otbConvolutionImageFilter.hxx.

◆ SetFilter()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual void otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::SetFilter ( ArrayType  filter)
inlinevirtual

Set the input filter

Definition at line 126 of file otbConvolutionImageFilter.h.

◆ SetNormalizeFilter()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual void otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::SetNormalizeFilter ( bool  _arg)
virtual

Set/Get methods for the normalization of the filter

◆ SetRadius()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
virtual void otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::SetRadius ( const InputSizeType  rad)
inlinevirtual

Set the radius of the neighborhood of the filter

Definition at line 104 of file otbConvolutionImageFilter.h.

◆ ThreadedGenerateData()

template<class TInputImage , class TOutputImage , class TBoundaryCondition , class TFilterPrecision >
void otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::ThreadedGenerateData ( const OutputImageRegionType outputRegionForThread,
itk::ThreadIdType  threadId 
)
overrideprotected

ConvolutionImageFilter can be implemented as a multithreaded filter. Therefore, this implementation provides a ThreadedGenerateData() routine which is called for each processing thread. The output image data is allocated automatically by the superclass prior to calling ThreadedGenerateData(). ThreadedGenerateData can only write to the portion of the output image specified by the parameter "outputRegionForThread"

See also
ImageToImageFilter::ThreadedGenerateData(), ImageToImageFilter::GenerateData()

Definition at line 95 of file otbConvolutionImageFilter.hxx.

◆ typedef()

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::typedef ( itk::Concept::HasNumericTraits< InputPixelType )

Begin concept checking This class requires InputHasNumericTraitsCheck in the form of (itk::Concept::HasNumericTraits<InputPixelType>)

Member Data Documentation

◆ InputImageDimension

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
const unsigned int otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::InputImageDimension = TInputImage::ImageDimension
static

Extract dimension from input and output image.

Definition at line 72 of file otbConvolutionImageFilter.h.

◆ m_Filter

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
ArrayType otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::m_Filter
private

Array containing the filter values

Definition at line 191 of file otbConvolutionImageFilter.h.

◆ m_NormalizeFilter

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
bool otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::m_NormalizeFilter
private

Flag for filter coefficients normalization

Definition at line 194 of file otbConvolutionImageFilter.h.

◆ m_Radius

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
InputSizeType otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::m_Radius
private

Radius of the filter

Definition at line 188 of file otbConvolutionImageFilter.h.

◆ OutputImageDimension

template<class TInputImage , class TOutputImage , class TBoundaryCondition = itk::ZeroFluxNeumannBoundaryCondition<TInputImage>, class TFilterPrecision = typename itk::NumericTraits<typename TInputImage::InternalPixelType>::RealType>
const unsigned int otb::ConvolutionImageFilter< TInputImage, TOutputImage, TBoundaryCondition, TFilterPrecision >::OutputImageDimension = TOutputImage::ImageDimension
static

Convenient typedefs for simplifying declarations.

Definition at line 73 of file otbConvolutionImageFilter.h.


The documentation for this class was generated from the following files:
otb::ConvolutionImageFilter::SetFilter
virtual void SetFilter(ArrayType filter)
Definition: otbConvolutionImageFilter.h:126