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
|
/*=========================================================================
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 otbBandMathImageFilter_h
#define otbBandMathImageFilter_h
#include "itkInPlaceImageFilter.h"
#include "itkImageRegionIteratorWithIndex.h"
#include "itkArray.h"
#include "otbParser.h"
namespace otb
{
/** \class BandMathImageFilter
* \brief Performs a mathematical operation on the input images
* according to the formula specified by the user.
*
* This filter is based on the mathematical parser library muParser.
* The built in functions and operators list is available at:
* http://muparser.sourceforge.net/mup_features.html#idDef2
*
* OTB additional functions:
* ndvi(r, niri)
*
* OTB additional constants:
* e - log2e - log10e - ln2 - ln10 - pi - euler
*
* In order to use this filter, at least one input image is to be
* set. An associated variable name can be specified or not by using
* the corresponding SetNthInput method. For the nth input image, if
* no associated variable name has been spefified, a default variable
* name is given by concatenating the letter "b" (for band) and the
* corresponding input index.
* Next step is to set the expression according to the variable
* names. For example, in the default case with three input images the
* following expression is valid :
* "ndvi(b1, b2)*b3"
*
* As an additional functionality, the filter also granted access to
* indexes information under special virtual bands named idxX, idxY
* for the images indexes and idxPhyX, idxPhyY for the physical
* indexes.
* It allows the user to perform, for example a spatial processing
* aiming to suppress a determined area :
* "if(sqrt((idxPhyX-105.3)*(idxPhyX-105.3)+
* (idxPhyY-207.1)*(idxPhyY-207.1))>100, b1, 0)"
* This expression replace the physical zone around the point of
* physical index (105.3; 207.1) by a black area
* This functionality assumes that all the band involved have the same
* spacing and origin.
*
*
* \sa Parser
*
* \ingroup Streamed
* \ingroup Threaded
*
* \ingroup OTBMathParser
*/
template< class TImage >
class ITK_EXPORT BandMathImageFilter
: public itk::InPlaceImageFilter< TImage >
{
public:
/** Standard class typedefs. */
typedef BandMathImageFilter< TImage > Self;
typedef itk::InPlaceImageFilter< TImage > Superclass;
typedef itk::SmartPointer< Self > Pointer;
typedef itk::SmartPointer< const Self > ConstPointer;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** Run-time type information (and related methods). */
itkTypeMacro(BandMathImageFilter, InPlaceImageFilter);
/** Some convenient typedefs. */
typedef TImage ImageType;
typedef typename ImageType::ConstPointer ImagePointer;
typedef typename ImageType::RegionType ImageRegionType;
typedef typename ImageType::PixelType PixelType;
typedef typename ImageType::IndexType IndexType;
typedef typename ImageType::PointType OrigineType;
typedef typename ImageType::SpacingType SpacingType;
typedef Parser ParserType;
typedef itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType;
/** Set the nth filter input with or without a specified associated variable name */
using Superclass::SetNthInput;
void SetNthInput( DataObjectPointerArraySizeType idx, const ImageType * image);
void SetNthInput( DataObjectPointerArraySizeType idx, const ImageType * image, const std::string& varName);
/** Change the nth filter input associated variable name */
void SetNthInputName(DataObjectPointerArraySizeType idx, const std::string& expression);
/** Set the expression to be parsed */
void SetExpression(const std::string& expression);
/** Return the expression to be parsed */
std::string GetExpression() const;
/** Return the nth filter input associated variable name */
std::string GetNthInputName(DataObjectPointerArraySizeType idx) const;
/** Return a pointer on the nth filter input */
ImageType * GetNthInput(DataObjectPointerArraySizeType idx);
protected :
BandMathImageFilter();
~BandMathImageFilter() ITK_OVERRIDE;
void PrintSelf(std::ostream& os, itk::Indent indent) const ITK_OVERRIDE;
void BeforeThreadedGenerateData() ITK_OVERRIDE;
void ThreadedGenerateData(const ImageRegionType& outputRegionForThread, itk::ThreadIdType threadId ) ITK_OVERRIDE;
void AfterThreadedGenerateData() ITK_OVERRIDE;
private :
BandMathImageFilter(const Self&); //purposely not implemented
void operator=(const Self&); //purposely not implemented
std::string m_Expression;
std::vector<ParserType::Pointer> m_VParser;
std::vector< std::vector<double> > m_AImage;
std::vector< std::string > m_VVarName;
unsigned int m_NbVar;
SpacingType m_Spacing;
OrigineType m_Origin;
long m_UnderflowCount;
long m_OverflowCount;
itk::Array<long> m_ThreadUnderflow;
itk::Array<long> m_ThreadOverflow;
};
}//end namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbBandMathImageFilter.txx"
#endif
#endif
|