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 "../srs/Converter.h"    43 #include "../geometry/LinearRing.h"    44 #include "../geometry/GTParameters.h"    45 #include "../geometry/MultiPoint.h"    46 #include "../geometry/Surface.h"    56 #include <boost/numeric/ublas/matrix.hpp>    57 #include <boost/shared_ptr.hpp>    84       SpectralSensorParams(
const int &band, 
const double &lower, 
const double &upper, 
const double &min, 
const double &max) :
   111       const std::vector< te::rst::BandProperty* >& bandsProperties,
   112       const std::string& outDataSetName,
   113       const std::string& dataSourceType,
   128       const std::vector< te::rst::BandProperty* >& bandsProperties,
   129       const std::string& outDataSetName,
   145       const std::vector< te::rst::BandProperty* >& bandsProperties,
   146       const std::map< std::string, std::string>& rasterInfo,
   147       const std::string& rasterType,
   148       std::unique_ptr< te::rst::Raster >& outRasterPtr );    
   160       std::vector< te::rst::BandProperty* > bandsProperties,
   174       std::vector< te::rst::BandProperty* > bandsProperties,
   175       const std::string& fileName,
   189       std::vector< te::rst::BandProperty* > bandsProperties,
   190       const std::string& fileName,
   191       std::unique_ptr< te::rst::Raster >& outRasterPtr );    
   205       const std::string& fileName );      
   230       unsigned int inputVectorSize, 
double* outputVector );
   242       unsigned int inputVectorSize, 
const int outputVectorDataType,
   243       void* outputVector );
   253     template< 
typename MatrixElementT >
   255       const bool normalize, 
const std::string& fileName )
   257       std::map<std::string, std::string> rInfo;
   258       rInfo[
"URI"] = fileName;
   260       std::vector<te::rst::BandProperty*> bandsProperties;
   266       bandsProperties[0]->m_noDataValue = -1.0 * std::numeric_limits<double>::max();
   271       std::unique_ptr< te::rst::Raster > outputRasterPtr(
   275       unsigned int line = 0;
   276       unsigned int col = 0;
   279       MatrixElementT matrixValue = 0;
   281       MatrixElementT gain = 1.0;
   282       MatrixElementT offset = 0.0;
   285         MatrixElementT matrixValueMin = std::numeric_limits< MatrixElementT >::max();
   286         MatrixElementT matrixValueMax = -1.0 * matrixValueMin;
   287         for( line = 0 ; line < nLines ; ++line )
   289           for( col = 0 ; col < nCols ; ++col )
   291             matrixValue = matrix[ line ][ col ];
   292             if( matrixValue < matrixValueMin )
   293               matrixValueMin = matrixValue;
   294             if( matrixValue > matrixValueMax )
   295               matrixValueMax = matrixValue;
   299         if( matrixValueMax == matrixValueMin )
   306           gain = 255.0 / ( matrixValueMax - matrixValueMin );
   307           offset = -1.0 * ( matrixValueMin );
   311       const MatrixElementT min0 = 0;
   312       const MatrixElementT max255 = 255;
   314       for( line = 0 ; line < nLines ; ++line )
   316         for( col = 0 ; col < nCols ; ++col )
   318           matrixValue = matrix[ line ][ col ];
   322             matrixValue += offset;
   324             matrixValue = std::max( min0, matrixValue );
   325             matrixValue = std::min( max255, matrixValue );
   328           outputRasterPtr->setValue( col, line, static_cast<double>(matrixValue), 0 );
   432       const unsigned int redBandIdx, 
const unsigned int greenBandIdx,
   433       const unsigned int blueBandIdx, 
const double rgbRangeMin, 
   453       const te::rst::Raster& inputGreenRaster, 
const unsigned int greenBandIdx,
   454       const te::rst::Raster& inputBlueRaster, 
const unsigned int blueBandIdx,
   455       const double rgbRangeMin, 
const double rgbRangeMax, 
te::rst::Raster& outputIHSRaster);
   471       const unsigned int intensityBandIdx, 
const unsigned int hueBandIdx,
   472       const unsigned int saturationBandIdx, 
const double rgbRangeMin, 
   492       const te::rst::Raster& inputSRaster, 
const unsigned int saturationBandIdx,
   493       const double rgbRangeMin, 
const double rgbRangeMax, 
te::rst::Raster& outputRGBRaster);
   509       const unsigned int redBandIdx, 
const unsigned int greenBandIdx,
   510       const unsigned int blueBandIdx, 
const double rgbRangeMin,
   529       const te::rst::Raster& inputGreenRaster, 
const unsigned int greenBandIdx,
   530       const te::rst::Raster& inputBlueRaster, 
const unsigned int blueBandIdx,
   531       const double rgbRangeMin, 
const double rgbRangeMax, 
te::rst::Raster& outputHLSRaster);
   546       const unsigned int hueBandIdx, 
const unsigned int lightBandIdx,
   547       const unsigned int saturationBandIdx, 
const double rgbRangeMin,
   566       const te::rst::Raster& inputSRaster, 
const unsigned int saturationBandIdx,
   567       const double rgbRangeMin, 
const double rgbRangeMax, 
te::rst::Raster& outputRGBRaster);
   581       const unsigned int inputBandIndex,
   582       const unsigned int maxThreads, 
   598       const unsigned int inputBandIndex,
   599       const unsigned int maxThreads, 
   600       double const * 
const meanValuePtr, 
   601       double& stdDevValue );         
   618       const unsigned int inputBandIndex1,
   620       const unsigned int inputBandIndex2,
   621       const bool useRasterCache,
   622       double const * 
const mean1ValuePtr, 
   623       double const * 
const mean2ValuePtr,  
   624       double& covarianceValue );      
   640       const std::vector< unsigned int >& inputRasterBands,
   641       boost::numeric::ublas::matrix< double >& pcaMatrix,
   643       const std::vector< unsigned int >& pcaRasterBands,
   644       const unsigned int maxThreads );  
   659       const boost::numeric::ublas::matrix< double >& pcaMatrix,
   661       const std::vector< unsigned int >& outputRasterBands,
   662       const unsigned int maxThreads );       
   679       const unsigned int& inputRasterBandIdx,
   680       const std::vector< std::pair< double, double > >& targetValues,
   681       const bool enableProgress,
   682       const std::vector< te::gm::Polygon* >& restrictionPols,
   684       const unsigned int& outputRasterBandIdx );
   700       const std::vector< unsigned int >& inputRasterBands,
   701       const boost::numeric::ublas::matrix< double >& remapMatrix,
   703       const std::vector< unsigned int >& outputRasterBands,
   704       const unsigned int maxThreads ); 
   718       const std::vector< unsigned int >& inputRasterBands,
   719       const std::vector< std::map<std::string, std::string> > & outputRastersInfos,
   720       const std::string& outputDataSourceType,
   721       std::vector< boost::shared_ptr< te::rst::Raster > > & outputRastersPtrs );   
   737       const std::vector< unsigned int >& inputRasterBands,
   739       const std::map<std::string, std::string>& outputRasterInfo,
   740       const std::string& outputDataSourceType,
   741       std::unique_ptr< te::rst::Raster >& outputRasterPtr );     
   769     TERPEXPORT boost::numeric::ublas::matrix< double > 
   786       const std::vector< unsigned int >& inputRasterBands,
   788       const unsigned int levelsNumber,
   789       const boost::numeric::ublas::matrix< double >& filter );     
   805       const unsigned int levelsNumber,
   807       const std::vector< unsigned int >& outputRasterBands ); 
   827       const std::vector< unsigned int >& inputRasterBands,
   829       const unsigned int firstRow,
   830       const unsigned int firstColumn, 
   831       const unsigned int height, 
   832       const unsigned int width,
   833       const unsigned int newheight, 
   834       const unsigned int newwidth, 
   835       const std::map<std::string, std::string>& rinfo,
   836       const std::string& dataSourceType,
   837       std::unique_ptr< te::rst::Raster >& resampledRasterPtr ); 
   845     template< 
typename ContainerT >
   847       const bool useTPSecondCoordPair )
   849       if( tiePoints.size() < 3 )
   857         typename ContainerT::const_iterator it =
   859         const typename ContainerT::const_iterator itE =
   864           if( useTPSecondCoordPair )
   872         std::unique_ptr< te::gm::Geometry > convexHullPolPtr( points.
convexHull() );
   874         if( dynamic_cast< te::gm::Surface* >( convexHullPolPtr.get() ) )
   891       const unsigned int paletteSize,
   892       const bool randomize,
   893       std::vector< te::rst::BandProperty::ColorEntry >& palette );
   913       const unsigned int inputRasterBand,
   914       const bool createPaletteRaster,
   915       const unsigned int slicesNumber,
   916       const bool eqHistogram,
   917       const std::map< std::string, std::string >& rasterInfo,
   918       const std::string& rasterType,
   919       const bool enableProgress,
   920       std::vector< te::rst::BandProperty::ColorEntry > 
const * 
const palettePtr,
   921       std::unique_ptr< te::rst::Raster >& outRasterPtr,
   922       std::map< double, double > limits);
   927 #endif  // __TERRALIB_RP_INTERNAL_FUNCTIONS_H 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...
 
Index into a lookup table. 
 
A raster band description. 
 
double GetTPConvexHullArea(const ContainerT &tiePoints, const bool useTPSecondCoordPair)
Returns the tie points converx hull area. 
 
TERPEXPORT void SaveSensorParams(std::map< std::string, SpectralSensorParams > &)
Saves in SpectralSensor.json file the spectral sensors parameters. 
 
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::map< double, double > limits)
Generate all wavelet planes from the given input raster. 
 
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. 
 
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)...
 
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 std::vector< std::string > GetBandNames()
Returns a vector os with band's names. 
 
TERASTEREXPORT void GetDataTypeRanges(const int &dataType, double &min, double &max)
Return the values range of a given data type. 
 
TERPEXPORT bool ComposeBands(te::rp::FeederConstRaster &feeder, const std::vector< unsigned int > &inputRasterBands, const te::rst::Interpolator::Method &interpMethod, const std::map< std::string, std::string > &outputRasterInfo, const std::string &outputDataSourceType, std::unique_ptr< te::rst::Raster > &outputRasterPtr)
Compose a set of bands into one multi-band raster. 
 
TERPEXPORT std::pair< double, double > GetDigitalNumberBandInfo(std::string bandName)
Returns the maximun and minimum digital numbers of a given sensor/band. 
 
An abstract class for data providers like a DBMS, Web Services or a regular file. ...
 
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 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 double GetDigitalNumberBandMin(std::string bandName)
Returns the minimum digital number of a given sensor/band. 
 
InterpolationMethod
Allowed interpolation methods. 
 
TERPEXPORT double GetSpectralBandMax(std::string bandName)
Returns the maximum reflectance value of a given sensor/band. 
 
TERPEXPORT bool GetMeanValue(const te::rst::Raster &inputRaster, const unsigned int inputBandIndex, const unsigned int maxThreads, double &meanValue)
Get the mean of band pixel values. 
 
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. 
 
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. 
 
#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...
 
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. 
 
bool CreateRasterFileFromMatrix(const te::rp::Matrix< MatrixElementT > &matrix, const bool normalize, const std::string &fileName)
Create a tiff file from a matrix. 
 
bool TERPEXPORT Copy2DiskRaster(const te::rst::Raster &inputRaster, const std::string &fileName)
Create a new raster into a GDAL datasource. 
 
A LinearRing is a LineString that is both closed and simple. 
 
MultiPoint is a GeometryCollection whose elements are restricted to points. 
 
TERPEXPORT std::map< std::string, SpectralSensorParams > getSensorParams()
Returns a map with spectral sensors parameters defined in SpectralSensor.json file. 
 
A point with x and y coordinate values. 
 
TERPEXPORT bool NormalizeRaster(te::rst::Raster &inputRaster, double nmin=0.0, double nmax=255.0)
Normalizes one raster in a given interval. 
 
An abstract class for raster data strucutures. 
 
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. 
 
WaveletAtrousFilterType
Wavelet Atrous Filter types. 
 
virtual Geometry * convexHull() const _NOEXCEPT_OP(false)
This method calculates the Convex Hull of a geometry. 
 
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. 
 
#define TERPEXPORT
You can use this macro in order to export/import classes and functions from this module. 
 
void TERPEXPORT Convert2DoublesVector(void *inputVector, const int inputVectorDataType, unsigned int inputVectorSize, double *outputVector)
Convert vector elements. 
 
TERPEXPORT bool GetDetailedExtent(const te::rst::Grid &grid, te::gm::LinearRing &detailedExtent)
Create a datailed extent from the given grid. 
 
SpectralSensorParams(const int &band, const double &lower, const double &upper, const double &min, const double &max)
 
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. 
 
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 boost::numeric::ublas::matrix< double > CreateWaveletAtrousFilter(const WaveletAtrousFilterType &filterType)
Create a Wavele Atrous Filter. 
 
unsigned int getLinesNumber() const
The number of current matrix lines. 
 
TERPEXPORT bool RasterResample(const te::rst::Raster &inputRaster, const std::vector< unsigned int > &inputRasterBands, const te::rst::Interpolator::Method interpMethod, 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. 
 
Configuration flags for the Raster Processing module of TerraLib. 
 
TERPEXPORT double GetSpectralBandMin(std::string bandName)
Returns the minimum reflectance value of a given sensor/band. 
 
TERPEXPORT bool GetCovarianceValue(const te::rst::Raster &inputRaster1, const unsigned int inputBandIndex1, const te::rst::Raster &inputRaster2, const unsigned int inputBandIndex2, const bool useRasterCache, double const *const mean1ValuePtr, double const *const mean2ValuePtr, double &covarianceValue)
Get the covariance of band pixel values. 
 
TERPEXPORT SpectralSensorParams GetSpectralBandInfo(std::string bandName)
Returns the maximun and minimum reflectance values of a given sensor/band. 
 
Feeder from a input rasters. 
 
void add(Geometry *g)
It adds the geometry into the collection. 
 
A generic template matrix. 
 
TERPEXPORT bool GetIndexedDetailedExtent(const te::rst::Grid &grid, te::gm::LinearRing &indexedDetailedExtent)
Create a indexed (lines,columns) datailed extent from the given grid. 
 
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 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 GetStdDevValue(const te::rst::Raster &inputRaster, const unsigned int inputBandIndex, const unsigned int maxThreads, double const *const meanValuePtr, double &stdDevValue)
Get the standard deviation of band pixel values. 
 
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. 
 
static Raster * make()
It creates and returns an empty raster with default raster driver. 
 
unsigned int getColumnsNumber() const
The number of current matrix columns. 
 
A rectified grid is the spatial support for raster data. 
 
TERPEXPORT std::string GetSensorFilename()
Returns a json filename with spectral sensors parameters. 
 
void TERPEXPORT ConvertDoublesVector(double *inputVector, unsigned int inputVectorSize, const int outputVectorDataType, void *outputVector)
Convert a doubles vector. 
 
Surface is an abstract class that represents a 2-dimensional geometric objects.