25#ifndef __TERRALIB_RP_INTERNAL_FUNCTIONS_H 
   26#define __TERRALIB_RP_INTERNAL_FUNCTIONS_H 
   34#include "../dataaccess/datasource/DataSource.h" 
   35#include "../raster/Raster.h" 
   36#include "../raster/RasterFactory.h" 
   37#include "../raster/Grid.h" 
   38#include "../raster/BandProperty.h" 
   39#include "../raster/RasterFactory.h" 
   40#include "../raster/Interpolator.h" 
   41#include "../raster/Utils.h" 
   42#include "../raster/PositionIterator.h" 
   43#include "../srs/Converter.h" 
   44#include "../geometry/LinearRing.h" 
   45#include "../geometry/GTParameters.h" 
   46#include "../geometry/MultiPoint.h" 
   47#include "../geometry/Surface.h" 
   57#include <boost/numeric/ublas/matrix.hpp> 
   58#include <boost/shared_ptr.hpp> 
   85      SpectralSensorParams(
const int &band, 
const double &lower, 
const double &upper, 
const double &min, 
const double &max) :
 
  112      const std::vector< te::rst::BandProperty* >& bandsProperties,
 
  113      const std::string& outDataSetName,
 
  114      const std::string& dataSourceType,
 
  129      const std::vector< te::rst::BandProperty* >& bandsProperties,
 
  130      const std::string& outDataSetName,
 
  146      const std::vector< te::rst::BandProperty* >& bandsProperties,
 
  147      const std::map< std::string, std::string>& rasterInfo,
 
  148      const std::string& rasterType,
 
  149      std::unique_ptr< te::rst::Raster >& outRasterPtr );    
 
  161      std::vector< te::rst::BandProperty* > bandsProperties,
 
  175      std::vector< te::rst::BandProperty* > bandsProperties,
 
  176      const std::string& fileName,
 
  190      std::vector< te::rst::BandProperty* > bandsProperties,
 
  191      const std::string& fileName,
 
  192      std::unique_ptr< te::rst::Raster >& outRasterPtr );    
 
  206      const std::string& fileName );      
 
  231      unsigned int inputVectorSize, 
double* outputVector );
 
  243      unsigned int inputVectorSize, 
const int outputVectorDataType,
 
  244      void* outputVector );
 
  254    template< 
typename MatrixElementT >
 
  256      const bool normalize, 
const std::string& fileName )
 
  258      std::map<std::string, std::string> rInfo;
 
  259      rInfo[
"URI"] = fileName;
 
  261      std::vector<te::rst::BandProperty*> bandsProperties;
 
  267      bandsProperties[0]->m_noDataValue = -1.0 * std::numeric_limits<double>::max();
 
  272      std::unique_ptr< te::rst::Raster > outputRasterPtr;
 
  275        outputRasterPtr.reset(
 
  280        outputRasterPtr.reset();
 
  284      unsigned int line = 0;
 
  285      unsigned int col = 0;
 
  288      MatrixElementT matrixValue = 0;
 
  290      MatrixElementT gain = 1.0;
 
  291      MatrixElementT offset = 0.0;
 
  294        MatrixElementT matrixValueMin = std::numeric_limits< MatrixElementT >::max();
 
  295        MatrixElementT matrixValueMax = -1.0 * matrixValueMin;
 
  296        for( line = 0 ; line < nLines ; ++line )
 
  298          for( col = 0 ; col < nCols ; ++col )
 
  300            matrixValue = matrix[ line ][ col ];
 
  301            if( matrixValue < matrixValueMin )
 
  302              matrixValueMin = matrixValue;
 
  303            if( matrixValue > matrixValueMax )
 
  304              matrixValueMax = matrixValue;
 
  308        if( matrixValueMax == matrixValueMin )
 
  315          gain = 255.0 / ( matrixValueMax - matrixValueMin );
 
  316          offset = -1.0 * ( matrixValueMin );
 
  320      const MatrixElementT min0 = 0;
 
  321      const MatrixElementT max255 = 255;
 
  323      for( line = 0 ; line < nLines ; ++line )
 
  325        for( col = 0 ; col < nCols ; ++col )
 
  327          matrixValue = matrix[ line ][ col ];
 
  331            matrixValue += offset;
 
  333            matrixValue = std::max( min0, matrixValue );
 
  334            matrixValue = std::min( max255, matrixValue );
 
  337          outputRasterPtr->setValue( col, line, 
static_cast<double>(matrixValue), 0 );
 
  383      std::map<std::string, struct SpectralSensorParams >& specBandInfos );
 
  452      const unsigned int redBandIdx, 
const unsigned int greenBandIdx,
 
  453      const unsigned int blueBandIdx, 
const double rgbRangeMin, 
 
  473      const te::rst::Raster& inputGreenRaster, 
const unsigned int greenBandIdx,
 
  474      const te::rst::Raster& inputBlueRaster, 
const unsigned int blueBandIdx,
 
  475      const double rgbRangeMin, 
const double rgbRangeMax, 
te::rst::Raster& outputIHSRaster);
 
  491      const unsigned int intensityBandIdx, 
const unsigned int hueBandIdx,
 
  492      const unsigned int saturationBandIdx, 
const double rgbRangeMin, 
 
  512      const te::rst::Raster& inputSRaster, 
const unsigned int saturationBandIdx,
 
  513      const double rgbRangeMin, 
const double rgbRangeMax, 
te::rst::Raster& outputRGBRaster);
 
  529      const unsigned int redBandIdx, 
const unsigned int greenBandIdx,
 
  530      const unsigned int blueBandIdx, 
const double rgbRangeMin,
 
  549      const te::rst::Raster& inputGreenRaster, 
const unsigned int greenBandIdx,
 
  550      const te::rst::Raster& inputBlueRaster, 
const unsigned int blueBandIdx,
 
  551      const double rgbRangeMin, 
const double rgbRangeMax, 
te::rst::Raster& outputHLSRaster);
 
  566      const unsigned int hueBandIdx, 
const unsigned int lightBandIdx,
 
  567      const unsigned int saturationBandIdx, 
const double rgbRangeMin,
 
  586      const te::rst::Raster& inputSRaster, 
const unsigned int saturationBandIdx,
 
  587      const double rgbRangeMin, 
const double rgbRangeMax, 
te::rst::Raster& outputRGBRaster);
 
  603      const unsigned int inputBandIndex,
 
  604      const unsigned int maxThreads, 
 
  605      const bool forceNoDataValue,
 
  606      const double noDataValue,
 
  624      const unsigned int inputBandIndex,
 
  625      const unsigned int maxThreads, 
 
  626      const bool forceNoDataValue,
 
  627      const double noDataValue,      
 
  628      double const * 
const meanValuePtr, 
 
  629      double& stdDevValue );         
 
  646      const unsigned int inputBandIndex1,
 
  648      const unsigned int inputBandIndex2,
 
  649      const unsigned int maxThreads,
 
  650      double const * 
const mean1ValuePtr, 
 
  651      double const * 
const mean2ValuePtr,
 
  652      double& covarianceValue );      
 
  668      const std::vector< unsigned int >& inputRasterBands,
 
  669      boost::numeric::ublas::matrix< double >& pcaMatrix,
 
  671      const std::vector< unsigned int >& pcaRasterBands,
 
  672      const unsigned int maxThreads );  
 
  687      const boost::numeric::ublas::matrix< double >& pcaMatrix,
 
  689      const std::vector< unsigned int >& outputRasterBands,
 
  690      const unsigned int maxThreads );       
 
  707      const unsigned int& inputRasterBandIdx,
 
  708      const std::vector< std::pair< double, double > >& targetValues,
 
  709      const bool enableProgress,
 
  710      const std::vector< te::gm::Polygon* >& restrictionPols,
 
  712      const unsigned int& outputRasterBandIdx );
 
  728      const std::vector< unsigned int >& inputRasterBands,
 
  729      const boost::numeric::ublas::matrix< double >& remapMatrix,
 
  731      const std::vector< unsigned int >& outputRasterBands,
 
  732      const unsigned int maxThreads ); 
 
  746      const std::vector< unsigned int >& inputRasterBands,
 
  747      const std::vector< std::map<std::string, std::string> > & outputRastersInfos,
 
  748      const std::string& outputDataSourceType,
 
  749      std::vector< boost::shared_ptr< te::rst::Raster > > & outputRastersPtrs );   
 
  767      const std::vector< te::rst::Raster const * >& inputRasters,
 
  768      const std::vector< unsigned int >& inputRasterBands,
 
  770      const std::map<std::string, std::string>& outputRasterInfo,
 
  771      const std::string& outputDataSourceType,
 
  772      const unsigned int maxThreads,
 
  773      const bool allowNoDataPixels,
 
  774      const bool enableProgress,
 
  775      std::unique_ptr< te::rst::Raster >& outputRasterPtr );
 
  803    TERPEXPORT boost::numeric::ublas::matrix< double > 
 
  820      const std::vector< unsigned int >& inputRasterBands,
 
  822      const unsigned int levelsNumber,
 
  823      const boost::numeric::ublas::matrix< double >& filter );     
 
  839      const unsigned int levelsNumber,
 
  841      const std::vector< unsigned int >& outputRasterBands ); 
 
  863      const std::vector< unsigned int >& inputRasterBands,
 
  865      const unsigned int interpWindowRadius, 
 
  866      const unsigned int firstRow,
 
  867      const unsigned int firstColumn, 
 
  868      const unsigned int height, 
 
  869      const unsigned int width,
 
  870      const unsigned int newheight, 
 
  871      const unsigned int newwidth, 
 
  872      const std::map<std::string, std::string>& rinfo,
 
  873      const std::string& dataSourceType,
 
  874      std::unique_ptr< te::rst::Raster >& resampledRasterPtr ); 
 
  883    template< 
typename ContainerT >
 
  885      const bool useTPSecondCoordPair )
 
  887      if( tiePoints.size() < 3 )
 
  895        typename ContainerT::const_iterator it =
 
  897        const typename ContainerT::const_iterator itE =
 
  902          if( useTPSecondCoordPair )
 
  910        std::unique_ptr< te::gm::Geometry > convexHullPolPtr( points.
convexHull() );
 
  930      const unsigned int paletteSize,
 
  931      const bool randomize,
 
  932      std::vector< te::rst::BandProperty::ColorEntry >& palette );
 
  953      const unsigned int inputRasterBand,
 
  954      const bool createPaletteRaster,
 
  955      const unsigned int slicesNumber,
 
  956      const bool eqHistogram,
 
  957      const std::map< std::string, std::string >& rasterInfo,
 
  958      const std::string& rasterType,
 
  959      const bool enableProgress,
 
  960      std::vector< te::rst::BandProperty::ColorEntry > 
const * 
const palettePtr,
 
  961      std::unique_ptr< te::rst::Raster >& outRasterPtr,
 
  962      std::vector< double > 
const* slicesLimitsPtr,
 
  963      std::vector< double > 
const* slicesOutputValuesPtr );
 
  977      const unsigned int& inputRasterBandIdx,
 
  978      const unsigned targetRow,
 
  979      const unsigned targetCol,
 
  980      const double& outputValue );
 
  996      const unsigned int& inputRasterBandIdx,
 
  998      const unsigned int& outputRasterBandIdx,
 
  999      const unsigned targetInputRow,
 
 1000      const unsigned targetInputCol,
 
 1001      const double& outputValue );
 
 1020      const unsigned int bandIndex,
 
 1021      const std::vector< te::gm::Polygon const* >& rois,
 
 1022      const unsigned int histoBins,
 
 1023      const unsigned int maxThreads,
 
 1025      std::map<double, unsigned>& rHistogram,
 
 1026      std::map<double, unsigned>& iHistogram );
 
 1047      const std::map<double, unsigned int>& histogram,
 
 1083      const unsigned int bandIndex1,
 
 1084      const unsigned int bandIndex2,      
 
 1085      const std::vector< te::gm::Polygon const* >& rois,
 
 1086      const unsigned int histoBins,
 
 1087      const unsigned int maxThreads,
 
 1089      std::map< std::pair< double, double >, 
unsigned int>& realJointHistogram,
 
 1090      std::map< std::pair< double, double >, 
unsigned int>& imagJointHistogram,
 
 1091      std::map<double, unsigned int> 
const * 
const rasterRealHistogram1ptr,
 
 1092      std::map<double, unsigned int> 
const * 
const rasterImagHistogram1Ptr,
 
 1093      std::map<double, unsigned int> 
const * 
const rasterRealHistogram2Ptr,
 
 1094      std::map<double, unsigned int> 
const * 
const rasterImagHistogram2Ptr );
 
 1107      const std::map< std::pair< double, double >, 
unsigned int>& jointHistogram,
 
 1108      const std::map<double, unsigned int>& histogram1,
 
 1109      const std::map<double, unsigned int>& histogram2,
 
 1111      double& corrlationCoef );
 
 1122      const unsigned int bandIndex,
 
 1123      const std::complex< double >& value );
 
 1139      const std::vector< unsigned int >& inputRasterBandsIdx,
 
 1140      const std::vector< double >& bandNoDataValues,
 
 1141      std::unique_ptr< te::gm::Polygon >& outPolygonPtr );    
 
#define TERP_TRUE_OR_RETURN_FALSE(value, message)
Checks if value is true. For false values a warning message will be logged and a return of context wi...
 
An abstract class for data providers like a DBMS, Web Services or a regular file.
 
void add(Geometry *g)
It adds the geometry into the collection.
 
virtual Geometry * convexHull() const _NOEXCEPT_OP(false)
This method calculates the Convex Hull of a geometry.
 
A LinearRing is a LineString that is both closed and simple.
 
MultiPoint is a GeometryCollection whose elements are restricted to points.
 
A point with x and y coordinate values.
 
Surface is an abstract class that represents a 2-dimensional geometric objects.
 
A generic template matrix.
 
unsigned int getLinesNumber() const
The number of current matrix lines.
 
unsigned int getColumnsNumber() const
The number of current matrix columns.
 
A raster band description.
 
A rectified grid is the spatial support for raster data.
 
static Raster * make()
It creates and returns an empty raster with default raster driver.
 
An abstract class for raster data strucutures.
 
TERPEXPORT bool RasterSlicing(const te::rst::Raster &inputRaster, const unsigned int inputRasterBand, const bool createPaletteRaster, const unsigned int slicesNumber, const bool eqHistogram, const std::map< std::string, std::string > &rasterInfo, const std::string &rasterType, const bool enableProgress, std::vector< te::rst::BandProperty::ColorEntry > const *const palettePtr, std::unique_ptr< te::rst::Raster > &outRasterPtr, std::vector< double > const *slicesLimitsPtr, std::vector< double > const *slicesOutputValuesPtr)
Create a new raster grouping pixel values following the number of slices and/or slice limits.
 
double GetTPConvexHullArea(const ContainerT &tiePoints, const bool useTPSecondCoordPair)
Returns the tie points converx hull area.
 
TERPEXPORT std::string GetSensorFilename()
Returns a json filename with spectral sensors parameters.
 
TERPEXPORT bool ComposeBands(const std::vector< te::rst::Raster const * > &inputRasters, const std::vector< unsigned int > &inputRasterBands, const te::rst::Interpolator::Method &interpMethod, const std::map< std::string, std::string > &outputRasterInfo, const std::string &outputDataSourceType, const unsigned int maxThreads, const bool allowNoDataPixels, const bool enableProgress, std::unique_ptr< te::rst::Raster > &outputRasterPtr)
Compose a set of bands into one multi-band raster.
 
TERPEXPORT bool DirectWaveletAtrous(const te::rst::Raster &inputRaster, const std::vector< unsigned int > &inputRasterBands, te::rst::Raster &waveletRaster, const unsigned int levelsNumber, const boost::numeric::ublas::matrix< double > &filter)
Generate all wavelet planes from the given input raster.
 
TERPEXPORT bool GetExternalValidDataPolygon(const te::rst::Raster &inputRaster, const std::vector< unsigned int > &inputRasterBandsIdx, const std::vector< double > &bandNoDataValues, std::unique_ptr< te::gm::Polygon > &outPolygonPtr)
Compute the external raster valid data area polygon.
 
void TERPEXPORT ConvertDoublesVector(double *inputVector, unsigned int inputVectorSize, const int outputVectorDataType, void *outputVector)
Convert a doubles vector.
 
bool TERPEXPORT CreateNewMemRaster(const te::rst::Grid &rasterGrid, std::vector< te::rst::BandProperty * > bandsProperties, RasterHandler &outRasterHandler)
Create a new raster into a new memory datasource.
 
TERPEXPORT double GetSpectralBandMin(std::string bandName)
Returns the minimum reflectance value of a given sensor/band.
 
bool TERPEXPORT CreateNewRaster(const te::rst::Grid &rasterGrid, const std::vector< te::rst::BandProperty * > &bandsProperties, const std::string &outDataSetName, const std::string &dataSourceType, RasterHandler &outRasterHandler)
Create a new raster into the givem data source.
 
TERPEXPORT std::map< std::string, SpectralSensorParams > getSensorParams()
Returns a map with spectral sensors parameters defined in SpectralSensor.json file.
 
TERPEXPORT bool GetStdDevValue(const te::rst::Raster &inputRaster, const unsigned int inputBandIndex, const unsigned int maxThreads, const bool forceNoDataValue, const double noDataValue, double const *const meanValuePtr, double &stdDevValue)
Get the standard deviation of band pixel values.
 
TERPEXPORT void SaveSensorParams(std::map< std::string, SpectralSensorParams > &)
Saves in SpectralSensor.json file the spectral sensors parameters.
 
TERPEXPORT boost::numeric::ublas::matrix< double > CreateWaveletAtrousFilter(const WaveletAtrousFilterType &filterType)
Create a Wavele Atrous Filter.
 
TERPEXPORT void getHistogramStats(const std::map< double, unsigned int > &histogram, double &min, double &max, double &mean, double &stdDev, double &mode, double &entropy, double &sum, double &sum2, double &sum3, double &sum4, double &variance, double &median)
Compute statiscts from the given histogram.
 
TERPEXPORT double GetSpectralBandMax(std::string bandName)
Returns the maximum reflectance value of a given sensor/band.
 
TERPEXPORT double GetDigitalNumberBandMin(std::string bandName)
Returns the minimum digital number of a given sensor/band.
 
bool TERPEXPORT Copy2DiskRaster(const te::rst::Raster &inputRaster, const std::string &fileName)
Create a new raster into a GDAL datasource.
 
TERPEXPORT bool ConvertIHS2RGB(const te::rst::Raster &inputIHSRaster, const unsigned int intensityBandIdx, const unsigned int hueBandIdx, const unsigned int saturationBandIdx, const double rgbRangeMin, const double rgbRangeMax, te::rst::Raster &outputRGBRaster)
IHS to RGB conversion.
 
TERPEXPORT bool InversePrincipalComponents(const te::rst::Raster &pcaRaster, const boost::numeric::ublas::matrix< double > &pcaMatrix, te::rst::Raster &outputRaster, const std::vector< unsigned int > &outputRasterBands, const unsigned int maxThreads)
Regenerate the original raster from its principal components.
 
TERPEXPORT bool ConvertRGB2HLS(const te::rst::Raster &inputRGBRaster, const unsigned int redBandIdx, const unsigned int greenBandIdx, const unsigned int blueBandIdx, const double rgbRangeMin, const double rgbRangeMax, te::rst::Raster &outputHLSRaster)
RGB to HLS conversion.
 
TERPEXPORT bool GetCovarianceValue(const te::rst::Raster &inputRaster1, const unsigned int inputBandIndex1, const te::rst::Raster &inputRaster2, const unsigned int inputBandIndex2, const unsigned int maxThreads, double const *const mean1ValuePtr, double const *const mean2ValuePtr, double &covarianceValue)
Get the covariance of band pixel values.
 
TERPEXPORT std::vector< std::string > GetBandNames()
Returns a vector os with band's names.
 
TERPEXPORT bool NormalizeRaster(te::rst::Raster &inputRaster, double nmin=0.0, double nmax=255.0)
Normalizes one raster in a given interval.
 
TERPEXPORT bool ReplaceContiguousSegmentValues(te::rst::Raster &inputRaster, const unsigned int &inputRasterBandIdx, const unsigned targetRow, const unsigned targetCol, const double &outputValue)
Replace a contiguos segment pixel values.
 
bool CreateRasterFileFromMatrix(const te::rp::Matrix< MatrixElementT > &matrix, const bool normalize, const std::string &fileName)
Create a tiff file from a matrix.
 
TERPEXPORT double GetDigitalNumberBandMax(std::string bandName)
Returns the maximum digital number of a given sensor/band.
 
TERPEXPORT bool SelectiveReplaceValues(const te::rst::Raster &inputRaster, const unsigned int &inputRasterBandIdx, const std::vector< std::pair< double, double > > &targetValues, const bool enableProgress, const std::vector< te::gm::Polygon * > &restrictionPols, te::rst::Raster &outputRaster, const unsigned int &outputRasterBandIdx)
Remap all pixel values using a user supplied target values (non-target values are just copied from in...
 
TERPEXPORT void GetAllSpectralBandInfos(std::map< std::string, struct SpectralSensorParams > &specBandInfos)
Get all spectral band infos.
 
TERPEXPORT bool getHistograms(const te::rst::Raster &raster, const unsigned int bandIndex, const std::vector< te::gm::Polygon const * > &rois, const unsigned int histoBins, const unsigned int maxThreads, const te::rst::PolygonIterator< double >::IterationType roisItType, std::map< double, unsigned > &rHistogram, std::map< double, unsigned > &iHistogram)
Compute and return the histogram soccurring values (real and imaginary) inside defined regions of int...
 
bool TERPEXPORT CreateNewGdalRaster(const te::rst::Grid &rasterGrid, std::vector< te::rst::BandProperty * > bandsProperties, const std::string &fileName, RasterHandler &outRasterHandler)
Create a new raster into a GDAL datasource.
 
TERPEXPORT bool GetMeanValue(const te::rst::Raster &inputRaster, const unsigned int inputBandIndex, const unsigned int maxThreads, const bool forceNoDataValue, const double noDataValue, double &meanValue)
Get the mean of band pixel values.
 
TERPEXPORT bool DirectPrincipalComponents(const te::rst::Raster &inputRaster, const std::vector< unsigned int > &inputRasterBands, boost::numeric::ublas::matrix< double > &pcaMatrix, te::rst::Raster &pcaRaster, const std::vector< unsigned int > &pcaRasterBands, const unsigned int maxThreads)
Generate all principal components from the given input raster.
 
TERPEXPORT bool FillBand(te::rst::Raster &raster, const unsigned int bandIndex, const std::complex< double > &value)
 
TERPEXPORT bool GetDetailedExtent(const te::rst::Grid &grid, te::gm::LinearRing &detailedExtent)
Create a datailed extent from the given grid.
 
void TERPEXPORT GetDataTypeRange(const int dataType, double &min, double &max)
Returns the real data type range (all values that can be represented by the given data type).
 
TERPEXPORT bool InverseWaveletAtrous(const te::rst::Raster &waveletRaster, const unsigned int levelsNumber, te::rst::Raster &outputRaster, const std::vector< unsigned int > &outputRasterBands)
Regenerate the original raster from its wavelets planes.
 
TERPEXPORT void CreateFixedStepPalette(const unsigned int paletteSize, const bool randomize, std::vector< te::rst::BandProperty::ColorEntry > &palette)
Create a fixed step sequential color palette.
 
TERPEXPORT bool ConvertRGB2IHS(const te::rst::Raster &inputRGBRaster, const unsigned int redBandIdx, const unsigned int greenBandIdx, const unsigned int blueBandIdx, const double rgbRangeMin, const double rgbRangeMax, te::rst::Raster &outputIHSRaster)
RGB to IHS conversion.
 
TERPEXPORT bool ConvertHLS2RGB(const te::rst::Raster &inputHLSRaster, const unsigned int hueBandIdx, const unsigned int lightBandIdx, const unsigned int saturationBandIdx, const double rgbRangeMin, const double rgbRangeMax, te::rst::Raster &outputRGBRaster)
HLS to RGB conversion.
 
TERPEXPORT bool RasterResample(const te::rst::Raster &inputRaster, const std::vector< unsigned int > &inputRasterBands, const te::rst::Interpolator::Method interpMethod, const unsigned int interpWindowRadius, const unsigned int firstRow, const unsigned int firstColumn, const unsigned int height, const unsigned int width, const unsigned int newheight, const unsigned int newwidth, const std::map< std::string, std::string > &rinfo, const std::string &dataSourceType, std::unique_ptr< te::rst::Raster > &resampledRasterPtr)
Resample a subset of the raster, given a box.
 
TERPEXPORT bool DecomposeBands(const te::rst::Raster &inputRaster, const std::vector< unsigned int > &inputRasterBands, const std::vector< std::map< std::string, std::string > > &outputRastersInfos, const std::string &outputDataSourceType, std::vector< boost::shared_ptr< te::rst::Raster > > &outputRastersPtrs)
Decompose a multi-band raster into a set of one-band rasters.
 
TERPEXPORT bool RemapValues(const te::rst::Raster &inputRaster, const std::vector< unsigned int > &inputRasterBands, const boost::numeric::ublas::matrix< double > &remapMatrix, te::rst::Raster &outputRaster, const std::vector< unsigned int > &outputRasterBands, const unsigned int maxThreads)
Remap pixel values using a remap function matrix.
 
WaveletAtrousFilterType
Wavelet Atrous Filter types.
 
TERPEXPORT bool getJointHistograms(const te::rst::Raster &raster, const unsigned int bandIndex1, const unsigned int bandIndex2, const std::vector< te::gm::Polygon const * > &rois, const unsigned int histoBins, const unsigned int maxThreads, const te::rst::PolygonIterator< double >::IterationType roisItType, std::map< std::pair< double, double >, unsigned int > &realJointHistogram, std::map< std::pair< double, double >, unsigned int > &imagJointHistogram, std::map< double, unsigned int > const *const rasterRealHistogram1ptr, std::map< double, unsigned int > const *const rasterImagHistogram1Ptr, std::map< double, unsigned int > const *const rasterRealHistogram2Ptr, std::map< double, unsigned int > const *const rasterImagHistogram2Ptr)
Compute and return the joint histogram of occurring values (real and imaginary) inside defined region...
 
TERPEXPORT bool GetSpectralBandInfo(const std::string &bandName, SpectralSensorParams &specBandInfo)
Get the maximun and minimum reflectance values of a given sensor/band.
 
TERPEXPORT bool GetIndexedDetailedExtent(const te::rst::Grid &grid, te::gm::LinearRing &indexedDetailedExtent)
Create a indexed (lines,columns) datailed extent from the given grid.
 
void TERPEXPORT Convert2DoublesVector(void *inputVector, const int inputVectorDataType, unsigned int inputVectorSize, double *outputVector)
Convert vector elements.
 
TERPEXPORT std::pair< double, double > GetDigitalNumberBandInfo(std::string bandName)
Returns the maximun and minimum digital numbers of a given sensor/band.
 
@ B3SplineWAFilter
Spline filter type.
 
@ TriangleWAFilter
Triangle filter type.
 
@ InvalidWAFilter
Invalid filter type.
 
TERPEXPORT bool getJointHistogramStats(const std::map< std::pair< double, double >, unsigned int > &jointHistogram, const std::map< double, unsigned int > &histogram1, const std::map< double, unsigned int > &histogram2, double &covariance, double &corrlationCoef)
Compute statiscts from the given histogram.
 
TERASTEREXPORT void GetDataTypeRanges(const int &dataType, double &min, double &max)
Return the values range of a given data type.
 
InterpolationMethod
Allowed interpolation methods.
 
@ GrayIdxCInt
Index into a lookup table.
 
SpectralSensorParams(const int &band, const double &lower, const double &upper, const double &min, const double &max)
 
#define TERPEXPORT
You can use this macro in order to export/import classes and functions from this module.
 
Proxy configuration file for TerraView (see terraview_config.h).