23 #ifndef otbBandMathXImageFilter_h 24 #define otbBandMathXImageFilter_h 26 #include "itkConstNeighborhoodIterator.h" 27 #include "itkImageToImageFilter.h" 64 template <
class TImage>
70 typedef itk::ImageToImageFilter<TImage, TImage>
Superclass;
85 typedef typename itk::ConstNeighborhoodIterator<TImage>::RadiusType
RadiusType;
101 using Superclass::SetNthInput;
102 void SetNthInput(DataObjectPointerArraySizeType idx,
const ImageType* image);
103 void SetNthInput(DataObjectPointerArraySizeType idx,
const ImageType* image,
const std::string& varName);
107 ImageType* GetNthInput(DataObjectPointerArraySizeType idx);
110 void SetManyExpressions(
bool flag);
113 void SetExpression(
const std::string& expression);
116 std::string GetExpression(
unsigned int IDExpression)
const;
119 void SetMatrix(
const std::string&
name,
const std::string& definition);
122 void SetConstant(
const std::string&
name,
double value);
125 void ExportContext(
const std::string& filename);
128 void ImportContext(
const std::string& filename);
131 void ClearExpression();
134 std::vector<std::string> GetVarNames()
const;
138 return !m_StatsVarDetected.empty();
144 void PrintSelf(std::ostream& os, itk::Indent indent)
const override;
146 void GenerateOutputInformation()
override;
147 void GenerateInputRequestedRegion()
override;
149 void BeforeThreadedGenerateData()
override;
150 void ThreadedGenerateData(
const ImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
override;
151 void AfterThreadedGenerateData()
override;
164 void operator=(
const Self&) =
delete;
167 void CheckImageDimensions();
168 void PrepareParsers();
169 void PrepareParsersGlobStats();
170 void OutputsDimensions();
173 std::vector<std::vector<ParserType::Pointer>>
m_VParser;
199 #ifndef OTB_MANUAL_INSTANTIATION std::vector< std::vector< adhocStruct > > m_AImage
bool GlobalStatsDetected() const
StreamingStatisticsVectorImageFilterType::MatrixType MatrixType
itk::Index< Monteverdi_DIMENSION > IndexType
static constexpr char const * name() noexcept
ImageType::ConstPointer ConstImagePointer
BandMathXImageFilter< TImage > Self
itk::Array< long > m_ThreadOverflow
ParserType::ValueType ValueType
itk::SmartPointer< Self > Pointer
std::vector< std::string > m_Expression
Definition of the standard floating point parser. Standard implementation of the mathematical express...
itk::SmartPointer< const Self > ConstPointer
StreamingStatisticsVectorImageFilterType::Pointer StreamingStatisticsVectorImageFilterPointerType
itk::SmartPointer< Self > Pointer
ImageType::PixelType PixelType
std::vector< unsigned int > m_outputsDimensions
itk::ConstNeighborhoodIterator< TImage >::RadiusType RadiusType
std::vector< adhocStruct > m_VNotAllowedVarName
ImageType::RegionType ImageRegionType
itk::ImageToImageFilter< TImage, TImage > Superclass
ImageType::PixelType::ValueType PixelValueType
std::vector< adhocStruct > m_VAllowedVarNameAuto
ImageType::SpacingType SpacingType
itk::Array< long > m_ThreadUnderflow
std::vector< adhocStruct > m_VAllowedVarNameAddedByUser
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
unsigned int m_SizeNeighbourhood
static constexpr GLenum value() noexcept
std::vector< int > m_StatsVarDetected
std::vector< unsigned int > m_NeighDetected
ImageType::Pointer ImagePointer
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
StreamingStatisticsVectorImageFilter< ImageType > StreamingStatisticsVectorImageFilterType
std::vector< std::vector< ParserType::Pointer > > m_VParser
std::vector< RadiusType > m_NeighExtremaSizes
std::vector< adhocStruct > m_VVarName
This class streams the whole input image through the PersistentStatisticsImageFilter.
VectorImageType::SpacingType SpacingType
ImageType::PointType OrigineType
std::vector< adhocStruct > m_VFinalAllowedVarName
StatFilterType::MatrixType MatrixType
VectorImageType::PointType PointType
ImageType::IndexType IndexType
Performs mathematical operations on the input images according to the formula specified by the user...