58#ifndef LIBRTCORE_H_INCLUDED
59#define LIBRTCORE_H_INCLUDED
71#if defined(__FreeBSD__) || defined(__FreeBSD_kernel__) || defined(__OpenBSD__)
95#if defined(__WIN32__) || defined(__NT__) || defined(_WIN32)
99#if defined(__BORLANDC__) && defined(MSDOS)
104#if defined(__APPLE__)
110#if defined(sun) || defined(__sun)
127#include "gdal_frmts.h"
128#include "gdalwarper.h"
131#include "cpl_string.h"
132#include "cpl_minixml.h"
134#include "ogr_srs_api.h"
136#include "postgis_config.h"
140# define __attribute__ (x)
233typedef char* (*rt_options)(
const char* varname);
234typedef void* (*rt_allocator)(
size_t size);
235typedef void* (*rt_reallocator)(
void *mem,
size_t size);
253extern void*
rtalloc(
size_t size);
254extern void*
rtrealloc(
void* mem,
size_t size);
285#if POSTGIS_DEBUG_LEVEL > 0
288#define RASTER_DEBUG(level, msg) \
290 if (POSTGIS_DEBUG_LEVEL >= level) \
291 rtinfo("[%s:%s:%d] " msg, __FILE__, __func__, __LINE__); \
295#define RASTER_DEBUGF(level, msg, ...) \
297 if (POSTGIS_DEBUG_LEVEL >= level) \
298 rtinfo("[%s:%s:%d] " msg, __FILE__, __func__, __LINE__, __VA_ARGS__); \
304#define RASTER_DEBUG(level, msg) \
308#define RASTER_DEBUGF(level, msg, ...) \
375 double val,
double refval,
406 uint16_t distancex, uint16_t distancey,
435 uint16_t width, uint16_t height,
437 uint32_t hasnodata,
double nodataval,
473 uint16_t width, uint16_t height,
475 uint32_t hasnodata,
double nodataval,
476 uint8_t bandNum,
const char* path
719 void *vals, uint32_t len
764 void **vals, uint16_t *nvals
806 uint16_t distancex, uint16_t distancey,
807 int exclude_nodata_value,
823 rt_band band,
int exclude_nodata_value,
824 double *searchset,
int searchcount,
874 double *newval,
int *corrected
892 int exclude_nodata_value,
double sample,
int inc_vals,
893 uint64_t *cK,
double *cM,
double *cQ
916 uint32_t bin_count,
double *bin_widths, uint32_t bin_widths_count,
917 int right,
double min,
double max,
933 double *quantiles,
int quantiles_count, uint32_t *rtn_count);
970 int exclude_nodata_value,
double sample,
973 double *quantiles, uint32_t quantiles_count,
992 rt_band band,
int exclude_nodata_value,
993 double *search_values, uint32_t search_values_count,
double roundto,
994 uint32_t *rtn_total, uint32_t *rtn_count
1011 uint32_t hasnodata,
double nodataval,
1026 uint32_t hasnodata,
double nodataval
1149 double initialvalue,
1150 uint32_t hasnodata,
double nodatavalue,
1164 double scaleX,
double scaleY);
1194 double x,
double y);
1224 double skewX,
double skewY);
1255 double *i_mag,
double *j_mag,
double *theta_i,
double *theta_ij) ;
1273 double i_mag,
double j_mag,
double theta_i,
double theta_ij) ;
1292 double xskew,
double yskew,
double yscale,
1293 double *i_mag,
double *j_mag,
double *theta_i,
double *theta_ij) ;
1313 double j_mag,
double theta_i,
double theta_ij,
1314 double *xscale,
double *xskew,
double *yskew,
double *yscale) ;
1343 double *gt,
double *igt
1379 double xr,
double yr,
1380 double* xw,
double* yw,
1398 double xw,
double yw,
1399 double *xr,
double *yr,
1418 double xw,
double yw,
1419 double *xr,
double *yr,
1474 double xr,
double yr,
1475 double *r_value,
int *r_nodata
1498 double xr,
double yr,
1500 double *r_value,
int *r_nodata
1626 int exclude_nodata_value,
1683 int fromindex,
int toindex
1735 char *format,
char **options, uint64_t *gdalsize);
1768 int *excludeNodataValues,
1770 GDALDriverH *rtn_drv,
int *destroy_rtn_drv
1796 const char* src_srs,
1797 double contour_interval,
1798 double contour_base,
1799 int fixed_level_count,
1800 double *fixed_levels,
1848 const char *src_srs,
const char *dst_srs,
1849 double *scale_x,
double *scale_y,
1850 int *width,
int *height,
1851 double *ul_xw,
double *ul_yw,
1852 double *grid_xw,
double *grid_yw,
1853 double *skew_x,
double *skew_y,
1854 GDALResampleAlg resample_alg,
double max_err);
1883 uint32_t wkb_len,
const char *srs,
1885 double *init,
double *value,
1886 double *nodata, uint8_t *hasnodata,
1887 int *width,
int *height,
1888 double *scale_x,
double *scale_y,
1889 double *ul_xw,
double *ul_yw,
1890 double *grid_xw,
double *grid_yw,
1891 double *skew_x,
double *skew_y,
2113 int *aligned,
char **reason
2178 uint8_t hasnodata,
double nodataval,
2179 uint16_t distancex, uint16_t distancey,
2207#if POSTGIS_GEOS_VERSION >= 31400
2220rt_raster rt_raster_intersection_fractions(
2231extern void *
rtalloc(
size_t size);
2232extern void *
rtrealloc(
void *mem,
size_t size);
2239#define GDAL_ENABLE_ALL "ENABLE_ALL"
2240#define GDAL_DISABLE_ALL "DISABLE_ALL"
2241#define GDAL_VSICURL "VSICURL"
2247#if !defined(POSTGIS_RASTER_WARN_ON_TRUNCATION)
2248#define POSTGIS_RASTER_WARN_ON_TRUNCATION 0
2251#define POSTGIS_RT_1BBMAX 1
2252#define POSTGIS_RT_2BUIMAX 3
2253#define POSTGIS_RT_4BUIMAX 15
2254#define POSTGIS_RT_16F_MAX 65504.0F
2294 double initialvalue,
2295 int32_t checkvalint, uint32_t checkvaluint,
2296 float checkvalfloat,
double checkvaldouble,
2435#define FLT_NEQ(x, y) ((x != y) && !(isnan(x) && isnan(y)) && (fabs(x - y) > FLT_EPSILON))
2436#define FLT_EQ(x, y) ((x == y) || (isnan(x) && isnan(y)) || (fabs(x - y) <= FLT_EPSILON))
2437#define DBL_NEQ(x, y) ((x != y) && !(isnan(x) && isnan(y)) && (fabs(x - y) > DBL_EPSILON))
2438#define DBL_EQ(x, y) ((x == y) || (isnan(x) && isnan(y)) || (fabs(x - y) <= DBL_EPSILON))
2443#define ROUND(x, y) (((x > 0.0) ? floor((x * pow(10, y) + 0.5)) : ceil((x * pow(10, y) - 0.5))) / pow(10, y))
void rt_band_init_value(rt_band band, double initval)
Fill in the cells of a band with a starting value frequently used to init with nodata value.
rt_errorstate rt_raster_contains(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *contains)
Return ES_ERROR if error occurred in function.
rt_band rt_band_new_inline(uint16_t width, uint16_t height, rt_pixtype pixtype, uint32_t hasnodata, double nodataval, uint8_t *data)
Create an in-db rt_band with no data.
void rt_band_set_ownsdata_flag(rt_band band, int flag)
void(* rt_message_handler)(const char *string, va_list ap) __attribute__((format(printf
rt_band rt_band_reclass(rt_band srcband, rt_pixtype pixtype, uint32_t hasnodata, double nodataval, rt_reclassexpr *exprset, int exprcount)
Returns new band with values reclassified.
int rt_band_get_pixel_of_value(rt_band band, int exclude_nodata_value, double *searchset, int searchcount, rt_pixel *pixels)
Search band for pixel(s) with search values.
void rterror(const char *fmt,...) __attribute__((format(printf
Wrappers used for reporting errors and info.
void * rtalloc(size_t size)
Wrappers used for managing memory.
rt_pixtype rt_util_gdal_datatype_to_pixtype(GDALDataType gdt)
Convert GDALDataType to rt_pixtype.
rt_band rt_band_duplicate(rt_band band)
Create a new band duplicated from source band.
uint16_t rt_band_get_width(rt_band band)
Return width of this band.
rt_errorstate rt_raster_fully_within_distance(rt_raster rast1, int nband1, rt_raster rast2, int nband2, double distance, int *dfwithin)
Return ES_ERROR if error occurred in function.
void rt_band_set_hasnodata_flag(rt_band band, int flag)
Set hasnodata flag value.
rt_errorstate rt_band_set_isnodata_flag(rt_band band, int flag)
Set isnodata flag value.
LWPOLY * rt_util_envelope_to_lwpoly(rt_envelope ext)
rt_errorstate rt_raster_cell_to_geopoint(rt_raster raster, double xr, double yr, double *xw, double *yw, double *gt)
Convert an xr, yr raster point to an xw, yw point on map.
int32_t rt_raster_get_srid(rt_raster raster)
Get raster's SRID.
void * default_rt_allocator(size_t size)
The default memory/logging handlers installed by lwgeom_install_default_allocators()
rt_errorstate rt_util_rgb_to_hsv(double rgb[3], double hsv[3])
int rt_util_gdal_register_all(int force_register_all)
double rt_raster_get_x_skew(rt_raster raster)
Get skew about the X axis.
double rt_raster_get_x_offset(rt_raster raster)
Get raster x offset, in projection units.
int rt_raster_calc_gt_coeff(double i_mag, double j_mag, double theta_i, double theta_ij, double *xscale, double *xskew, double *yskew, double *yscale)
Calculates the coefficients of a geotransform given the physically significant parameters describing ...
int rt_raster_generate_new_band(rt_raster raster, rt_pixtype pixtype, double initialvalue, uint32_t hasnodata, double nodatavalue, int index)
Generate a new inline band and add it to a raster.
int rt_util_gdal_configured(void)
rt_raster rt_raster_from_wkb(const uint8_t *wkb, uint32_t wkbsize)
Construct an rt_raster from a binary WKB representation.
LWPOLY * rt_raster_pixel_as_polygon(rt_raster raster, int x, int y)
Get a raster pixel as a polygon.
void rt_raster_set_geotransform_matrix(rt_raster raster, double *gt)
Set raster's geotransform using 6-element array.
int8_t rt_util_clamp_to_8BSI(double value)
struct rt_quantile_t * rt_quantile
uint8_t rt_util_clamp_to_1BB(double value)
void void rtinfo(const char *fmt,...) __attribute__((format(printf
float rt_util_clamp_to_16F(double value)
rt_raster rt_raster_colormap(rt_raster raster, int nband, rt_colormap colormap)
Returns a new raster with up to four 8BUI bands (RGBA) from applying a colormap to the user-specified...
rt_errorstate rt_band_load_offline_data(rt_band band)
Load offline band's data.
void void void rtwarn(const char *fmt,...) __attribute__((format(printf
rt_errorstate rt_raster_contains_properly(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *contains)
Return ES_ERROR if error occurred in function.
rt_gdaldriver rt_raster_gdal_drivers(uint32_t *drv_count, uint8_t cancc)
Returns a set of available GDAL drivers.
int32_t rt_util_clamp_to_32BSI(double value)
uint16_t rt_util_float_to_float16(float value)
struct rt_reclassexpr_t * rt_reclassexpr
void rt_raster_set_scale(rt_raster raster, double scaleX, double scaleY)
Set scale in projection units.
int rt_raster_add_band(rt_raster raster, rt_band band, int index)
Add band data to a raster.
rt_errorstate rt_raster_get_inverse_geotransform_matrix(rt_raster raster, double *gt, double *igt)
Get 6-element array of raster inverse geotransform matrix.
void(* rt_deallocator)(void *mem)
float rt_util_float16_to_float(uint16_t value)
int rt_band_get_hasnodata_flag(rt_band band)
Get hasnodata flag value.
rt_pixtype rt_pixtype_index_from_name(const char *pixname)
uint8_t * rt_raster_to_gdal(rt_raster raster, const char *srs, char *format, char **options, uint64_t *gdalsize)
Return formatted GDAL raster from raster.
void * rt_band_get_data(rt_band band)
Get pointer to raster band data.
rt_errorstate rt_raster_geopoint_to_cell(rt_raster raster, double xw, double yw, double *xr, double *yr, double *igt)
Convert an xw, yw map point to a xr, yr raster point.
rt_errorstate rt_band_get_pixel(rt_band band, int x, int y, double *value, int *nodata)
Get pixel value.
struct rt_mask_t * rt_mask
void rt_raster_destroy(rt_raster raster)
Release memory associated to a raster.
void rt_raster_set_skews(rt_raster raster, double skewX, double skewY)
Set skews about the X and Y axis.
int rt_band_get_isnodata_flag(rt_band band)
Get isnodata flag value.
rt_errorstate rt_raster_surface(rt_raster raster, int nband, LWMPOLY **surface)
Get a raster as a surface (multipolygon).
rt_errorstate rt_raster_geopoint_to_rasterpoint(rt_raster raster, double xw, double yw, double *xr, double *yr, double *igt)
Convert an xw,yw map point to a xr,yr raster point.
rt_raster rt_raster_new(uint32_t width, uint32_t height)
Construct a raster with given dimensions.
char * rt_util_gdal_convert_sr(const char *srs, int proj4)
rt_extenttype rt_util_extent_type(const char *name)
int rt_util_dbl_trunc_warning(double initialvalue, int32_t checkvalint, uint32_t checkvaluint, float checkvalfloat, double checkvaldouble, rt_pixtype pixtype)
GDALDataType rt_util_pixtype_to_gdal_datatype(rt_pixtype pt)
Convert rt_pixtype to GDALDataType.
rt_errorstate rt_raster_intersects(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *intersects)
Return ES_ERROR if error occurred in function.
rt_errorstate rt_raster_coveredby(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *coveredby)
Return ES_ERROR if error occurred in function.
double rt_pixtype_get_min_value(rt_pixtype pixtype)
Return minimum value possible for pixel type.
int rt_raster_has_band(rt_raster raster, int nband)
Return TRUE if the raster has a band of this number.
rt_raster rt_raster_compute_skewed_raster(rt_envelope extent, double *skew, double *scale, double tolerance)
struct rt_pixel_t * rt_pixel
void void void char * default_rt_options(const char *varname)
rt_errorstate rt_raster_touches(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *touches)
Return ES_ERROR if error occurred in function.
rt_raster rt_raster_from_gdal_dataset(GDALDatasetH ds)
Return a raster from a GDAL dataset.
int rt_util_gdal_progress_func(double dfComplete, const char *pszMessage, void *pProgressArg)
rt_quantile rt_band_get_quantiles_stream(rt_band band, int exclude_nodata_value, double sample, uint64_t cov_count, struct quantile_llist **qlls, uint32_t *qlls_count, double *quantiles, uint32_t quantiles_count, uint32_t *rtn_count)
Compute the default set of or requested quantiles for a coverage.
struct rt_colormap_entry_t * rt_colormap_entry
rt_geomval rt_raster_gdal_polygonize(rt_raster raster, int nband, int exclude_nodata_value, int *pnElements)
Returns a set of "geomval" value, one for each group of pixel sharing the same value for the provided...
double rt_raster_get_x_scale(rt_raster raster)
Get scale X in projection units.
rt_band rt_raster_replace_band(rt_raster raster, rt_band band, int index)
Replace band at provided index with new band.
uint8_t rt_util_clamp_to_2BUI(double value)
void default_rt_deallocator(void *mem)
const char * rt_pixtype_name(rt_pixtype pixtype)
int rt_pixtype_alignment(rt_pixtype pixtype)
Return alignment requirements for data in the given pixel type.
void * default_rt_reallocator(void *mem, size_t size)
rt_errorstate rt_pixtype_compare_clamped_values(rt_pixtype pixtype, double val, double refval, int *isequal)
Test to see if two values are equal when clamped.
const char * rt_util_gdal_version(const char *request)
uint8_t rt_util_clamp_to_8BUI(double value)
double rt_band_get_min_value(rt_band band)
Returns the minimal possible value for the band according to the pixel type.
uint64_t rt_band_get_file_size(rt_band band)
Return file size in bytes.
struct rt_bandstats_t * rt_bandstats
GDALDatasetH rt_util_gdal_open(const char *fn, GDALAccess fn_access, int shared)
int quantile_llist_destroy(struct quantile_llist **list, uint32_t list_count)
int rt_band_check_is_nodata(rt_band band)
Returns TRUE if the band is only nodata values.
void void void char * rtoptions(const char *varname)
Wrappers used for options.
rt_errorstate rt_band_set_nodata(rt_band band, double val, int *converted)
Set nodata value.
char * rt_raster_to_hexwkb(rt_raster raster, int outasin, uint32_t *hexwkbsize)
Return this raster in HEXWKB form (null-terminated hex)
rt_errorstate rt_band_set_pixel(rt_band band, int x, int y, double val, int *converted)
Set single pixel's value.
rt_errorstate
Enum definitions.
rt_errorstate rt_band_set_pixel_line(rt_band band, int x, int y, void *vals, uint32_t len)
Set values of multiple pixels.
rt_raster rt_raster_from_hexwkb(const char *hexwkb, uint32_t hexwkbsize)
Construct an rt_raster from a text HEXWKB representation.
rt_raster rt_raster_gdal_rasterize(const unsigned char *wkb, uint32_t wkb_len, const char *srs, uint32_t num_bands, rt_pixtype *pixtype, double *init, double *value, double *nodata, uint8_t *hasnodata, int *width, int *height, double *scale_x, double *scale_y, double *ul_xw, double *ul_yw, double *grid_xw, double *grid_yw, double *skew_x, double *skew_y, char **options)
Return a raster of the provided geometry.
void rt_band_destroy(rt_band band)
Destroy a raster band.
rt_errorstate rt_raster_overlaps(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *overlaps)
Return ES_ERROR if error occurred in function.
int rt_util_gdal_supported_sr(const char *srs)
void rt_raster_get_phys_params(rt_raster rast, double *i_mag, double *j_mag, double *theta_i, double *theta_ij)
Calculates and returns the physically significant descriptors embodied in the geotransform attached t...
int16_t rt_util_clamp_to_16BSI(double value)
rt_errorstate rt_band_corrected_clamped_value(rt_band band, double val, double *newval, int *corrected)
Correct value when clamped value is equal to clamped NODATA value.
void * rtrealloc(void *mem, size_t size)
GDALResampleAlg rt_util_gdal_resample_alg(const char *algname)
Convert cstring name to GDAL Resample Algorithm.
uint16_t rt_raster_get_num_bands(rt_raster raster)
struct rt_valuecount_t * rt_valuecount
struct rt_gdaldriver_t * rt_gdaldriver
rt_raster rt_raster_clone(rt_raster raster, uint8_t deep)
Clone an existing raster.
rt_errorstate rt_raster_get_convex_hull(rt_raster raster, LWGEOM **hull)
Get raster's convex hull.
uint16_t rt_raster_get_height(rt_raster raster)
struct rt_reclassmap_t * rt_reclassmap
void *(* rt_reallocator)(void *mem, size_t size)
int rt_util_gdal_driver_registered(const char *drv)
rt_band rt_band_new_offline(uint16_t width, uint16_t height, rt_pixtype pixtype, uint32_t hasnodata, double nodataval, uint8_t bandNum, const char *path)
Create an out-db rt_band.
rt_errorstate rt_band_get_ext_band_num(rt_band band, uint8_t *bandnum)
Return bands' external band number (only valid when rt_band_is_offline returns non-zero).
rt_errorstate rt_raster_iterator(rt_iterator itrset, uint16_t itrcount, rt_extenttype extenttype, rt_raster customextent, rt_pixtype pixtype, uint8_t hasnodata, double nodataval, uint16_t distancex, uint16_t distancey, rt_mask mask, void *userarg, int(*callback)(rt_iterator_arg arg, void *userarg, double *value, int *nodata), rt_raster *rtnraster)
n-raster iterator.
uint64_t rt_band_get_file_timestamp(rt_band band)
Return file timestamp.
struct rt_iterator_arg_t * rt_iterator_arg
void rt_raster_set_srid(rt_raster raster, int32_t srid)
Set raster's SRID.
rt_band rt_band_reclass_exact(rt_band srcband, rt_reclassmap map, uint32_t hasnodata, double nodataval)
Returns new band with values reclassified.
void rt_raster_set_phys_params(rt_raster rast, double i_mag, double j_mag, double theta_i, double theta_ij)
Calculates the geotransform coefficients and applies them to the supplied raster.
rt_errorstate rt_band_get_pixel_resample(rt_band band, double xr, double yr, rt_resample_type resample, double *r_value, int *r_nodata)
Retrieve a point value from the raster using a world coordinate and selected resampling method.
uint8_t * rt_raster_to_wkb(rt_raster raster, int outasin, uint32_t *wkbsize)
Return this raster in WKB form.
rt_raster rt_raster_from_band(rt_raster raster, uint32_t *bandNums, int count)
Construct a new rt_raster from an existing rt_raster and an array of band numbers.
struct rt_geomval_t * rt_geomval
rt_errorstate rt_band_get_nodata(rt_band band, double *nodata)
Get NODATA value.
uint8_t rt_util_clamp_to_4BUI(double value)
rt_errorstate rt_util_hsv_to_rgb(double hsv[3], double rgb[3])
void rt_util_to_ogr_envelope(rt_envelope ext, OGREnvelope *env)
rt_pixtype rt_band_get_pixtype(rt_band band)
Return pixeltype of this band.
char * rtstrdup(const char *str)
rt_errorstate rt_pixel_set_to_array(rt_pixel npixel, uint32_t count, rt_mask mask, int x, int y, uint16_t distancex, uint16_t distancey, double ***value, int ***nodata, int *dimx, int *dimy)
char *(* rt_options)(const char *varname)
Global functions for memory/logging handlers.
void *(* rt_allocator)(size_t size)
LWPOINT * rt_raster_pixel_as_centroid_point(rt_raster rast, int x, int y)
Get a raster pixel centroid point.
const char * rt_band_get_ext_path(rt_band band)
Return band's external path (only valid when rt_band_is_offline returns non-zero).
rt_band rt_band_new_offline_from_path(uint16_t width, uint16_t height, int hasnodata, double nodataval, uint8_t bandNum, const char *path, int force)
Create an out-db rt_band from path.
rt_histogram rt_band_get_histogram(rt_bandstats stats, uint32_t bin_count, double *bin_widths, uint32_t bin_widths_count, int right, double min, double max, uint32_t *rtn_count)
Count the distribution of data.
uint16_t rt_raster_get_width(rt_raster raster)
uint32_t rt_band_get_nearest_pixel(rt_band band, int x, int y, uint16_t distancex, uint16_t distancey, int exclude_nodata_value, rt_pixel *npixels)
Get nearest pixel(s) with value (not NODATA) to specified pixel.
void default_rt_error_handler(const char *fmt, va_list ap) __attribute__((format(printf
void rtdealloc(void *mem)
rt_errorstate rt_band_get_pixel_bilinear(rt_band band, double xr, double yr, double *r_value, int *r_nodata)
Retrieve a point value from the raster using a world coordinate and bilinear interpolation.
rt_errorstate rt_util_gdal_sr_auth_info(GDALDatasetH hds, char **authname, char **authcode)
Get auth name and code.
uint16_t rt_util_clamp_to_16BUI(double value)
rt_errorstate rt_raster_get_envelope_geom(rt_raster raster, LWGEOM **env)
Get raster's envelope as a geometry.
struct rt_colormap_t * rt_colormap
struct rt_iterator_t * rt_iterator
int rt_band_clamped_value_is_nodata(rt_band band, double val)
Compare clamped value to band's clamped NODATA value.
void * rt_raster_serialize(rt_raster raster)
Return this raster in serialized form.
uint32_t rt_util_clamp_to_32BUI(double value)
rt_raster rt_raster_gdal_warp(rt_raster raster, const char *src_srs, const char *dst_srs, double *scale_x, double *scale_y, int *width, int *height, double *ul_xw, double *ul_yw, double *grid_xw, double *grid_yw, double *skew_x, double *skew_y, GDALResampleAlg resample_alg, double max_err)
Return a warped raster using GDAL Warp API.
int rt_band_get_ownsdata_flag(rt_band band)
Return 0 (FALSE) or non-zero (TRUE) indicating if rt_band is responsible for managing the memory for ...
GDALDatasetH rt_raster_to_gdal_mem(rt_raster raster, const char *srs, uint32_t *bandNums, int *excludeNodataValues, int count, GDALDriverH *rtn_drv, int *destroy_rtn_drv)
Return GDAL dataset using GDAL MEM driver from raster.
void rt_set_handlers_options(rt_allocator allocator, rt_reallocator reallocator, rt_deallocator deallocator, rt_message_handler error_handler, rt_message_handler info_handler, rt_message_handler warning_handler, rt_options options_handler)
rt_errorstate rt_raster_get_perimeter(rt_raster raster, int nband, LWGEOM **perimeter)
Get raster perimeter.
int rt_band_is_offline(rt_band band)
Return non-zero if the given band data is on the filesystem.
rt_errorstate rt_raster_within_distance(rt_raster rast1, int nband1, rt_raster rast2, int nband2, double distance, int *dwithin)
Return ES_ERROR if error occurred in function.
rt_valuecount rt_band_get_value_count(rt_band band, int exclude_nodata_value, double *search_values, uint32_t search_values_count, double roundto, uint32_t *rtn_total, uint32_t *rtn_count)
Count the number of times provided value(s) occur in the band.
rt_errorstate rt_raster_from_two_rasters(rt_raster rast1, rt_raster rast2, rt_extenttype extenttype, rt_raster *rtnraster, double *offset)
void(*) voi rt_install_default_allocators)(void)
Apply the default memory management (malloc() and free()) and error handlers.
rt_geos_spatial_test
GEOS spatial relationship tests available.
struct rt_raster_t * rt_raster
Types definitions.
int rt_util_same_geotransform_matrix(double *gt1, double *gt2)
int rt_raster_gdal_contour(rt_raster src_raster, int src_band, int src_srid, const char *src_srs, double contour_interval, double contour_base, int fixed_level_count, double *fixed_levels, int polygonize, size_t *ncontours, struct rt_contour_t **contours)
Return palloc'ed list of contours.
int rt_raster_copy_band(rt_raster torast, rt_raster fromrast, int fromindex, int toindex)
Copy one band from one raster to another.
double rt_raster_get_y_scale(rt_raster raster)
Get scale Y in projection units.
float rt_util_clamp_to_32F(double value)
rt_errorstate rt_raster_get_envelope(rt_raster raster, rt_envelope *env)
Get raster's envelope.
rt_errorstate rt_raster_same_alignment(rt_raster rast1, rt_raster rast2, int *aligned, char **reason)
rt_errorstate rt_raster_covers(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *covers)
Return ES_ERROR if error occurred in function.
double rt_raster_get_y_skew(rt_raster raster)
Get skew about the Y axis.
void rt_raster_get_geotransform_matrix(rt_raster raster, double *gt)
Get 6-element array of raster geotransform matrix.
rt_quantile rt_band_get_quantiles(rt_bandstats stats, double *quantiles, int quantiles_count, uint32_t *rtn_count)
Compute the default set of or requested quantiles for a set of data the quantile formula used is same...
void void void default_rt_info_handler(const char *fmt, va_list ap) __attribute__((format(printf
rt_errorstate rt_raster_copy_to_geometry(rt_raster raster, uint32_t bandnum, char dim, rt_resample_type resample, const LWGEOM *lwgeom_in, LWGEOM **lwgeom_out)
Copy values from a raster to the points on a geometry using the requested interpolation type.
void rt_raster_set_offsets(rt_raster raster, double x, double y)
Set insertion points in projection units.
rt_bandstats rt_band_get_summary_stats(rt_band band, int exclude_nodata_value, double sample, int inc_vals, uint64_t *cK, double *cM, double *cQ)
Compute summary statistics for a band.
void rt_set_handlers(rt_allocator allocator, rt_reallocator reallocator, rt_deallocator deallocator, rt_message_handler error_handler, rt_message_handler info_handler, rt_message_handler warning_handler)
This function is called when the PostgreSQL backend is taking care of the memory and we want to use p...
void rt_util_from_ogr_envelope(OGREnvelope env, rt_envelope *ext)
struct rt_histogram_t * rt_histogram
uint16_t rt_band_get_height(rt_band band)
Return height of this band.
struct rt_band_t * rt_band
rt_raster rt_raster_deserialize(void *serialized, int header_only)
Return a raster from a serialized form.
int rt_pixtype_size(rt_pixtype pixtype)
Return size in bytes of a value in the given pixtype.
void rt_raster_calc_phys_params(double xscale, double xskew, double yskew, double yscale, double *i_mag, double *j_mag, double *theta_i, double *theta_ij)
Calculates the physically significant descriptors embodied in an arbitrary geotransform.
rt_errorstate rt_band_get_pixel_line(rt_band band, int x, int y, uint16_t len, void **vals, uint16_t *nvals)
Get values of multiple pixels.
double rt_raster_get_y_offset(rt_raster raster)
Get raster y offset, in projection units.
int rt_raster_is_empty(rt_raster raster)
Return TRUE if the raster is empty.
rt_band rt_raster_get_band(rt_raster raster, int bandNum)
Return Nth band, or NULL if unavailable.
void void default_rt_warning_handler(const char *fmt, va_list ap) __attribute__((format(printf
Datum contains(PG_FUNCTION_ARGS)
Datum coveredby(PG_FUNCTION_ARGS)
Datum covers(PG_FUNCTION_ARGS)
Datum touches(PG_FUNCTION_ARGS)
Datum overlaps(PG_FUNCTION_ARGS)
static double distance(double x1, double y1, double x2, double y2)
struct quantile_llist_element * prev
struct quantile_llist_element * next
struct quantile_llist_element * element
struct quantile_llist_element * head
struct quantile_llist_element * tail
struct quantile_llist_index * index
union rt_band_t::@12 data
enum rt_colormap_t::@13 method
struct rt_reclassexpr_t::rt_reclassrange src
struct rt_reclassexpr_t::rt_reclassrange dst
struct rt_classpair_t * pairs