GDAL
gdal_alg.h
Go to the documentation of this file.
1/******************************************************************************
2 * $Id: gdal_alg.h 8e5eeb35bf76390e3134a4ea7076dab7d478ea0e 2018-11-14 22:55:13 +0100 Even Rouault $
3 *
4 * Project: GDAL Image Processing Algorithms
5 * Purpose: Prototypes, and definitions for various GDAL based algorithms.
6 * Author: Frank Warmerdam, warmerdam@pobox.com
7 *
8 ******************************************************************************
9 * Copyright (c) 2001, Frank Warmerdam
10 * Copyright (c) 2008-2012, Even Rouault <even dot rouault at mines-paris dot org>
11 *
12 * Permission is hereby granted, free of charge, to any person obtaining a
13 * copy of this software and associated documentation files (the "Software"),
14 * to deal in the Software without restriction, including without limitation
15 * the rights to use, copy, modify, merge, publish, distribute, sublicense,
16 * and/or sell copies of the Software, and to permit persons to whom the
17 * Software is furnished to do so, subject to the following conditions:
18 *
19 * The above copyright notice and this permission notice shall be included
20 * in all copies or substantial portions of the Software.
21 *
22 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
23 * OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
24 * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
25 * THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
26 * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
27 * FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
28 * DEALINGS IN THE SOFTWARE.
29 ****************************************************************************/
30
31#ifndef GDAL_ALG_H_INCLUDED
32#define GDAL_ALG_H_INCLUDED
33
40#ifndef DOXYGEN_SKIP
41#include "gdal.h"
42#include "cpl_minixml.h"
43#include "ogr_api.h"
44#endif
45
47
48int CPL_DLL CPL_STDCALL GDALComputeMedianCutPCT( GDALRasterBandH hRed,
49 GDALRasterBandH hGreen,
50 GDALRasterBandH hBlue,
51 int (*pfnIncludePixel)(int,int,void*),
52 int nColors,
53 GDALColorTableH hColorTable,
54 GDALProgressFunc pfnProgress,
55 void * pProgressArg );
56
57int CPL_DLL CPL_STDCALL GDALDitherRGB2PCT( GDALRasterBandH hRed,
58 GDALRasterBandH hGreen,
59 GDALRasterBandH hBlue,
60 GDALRasterBandH hTarget,
61 GDALColorTableH hColorTable,
62 GDALProgressFunc pfnProgress,
63 void * pProgressArg );
64
65int CPL_DLL CPL_STDCALL GDALChecksumImage( GDALRasterBandH hBand,
66 int nXOff, int nYOff, int nXSize, int nYSize );
67
68CPLErr CPL_DLL CPL_STDCALL
70 GDALRasterBandH hProximityBand,
71 char **papszOptions,
72 GDALProgressFunc pfnProgress,
73 void * pProgressArg );
74
75CPLErr CPL_DLL CPL_STDCALL
77 GDALRasterBandH hMaskBand,
78 double dfMaxSearchDist,
79 int bDeprecatedOption,
80 int nSmoothingIterations,
81 char **papszOptions,
82 GDALProgressFunc pfnProgress,
83 void * pProgressArg );
84
85CPLErr CPL_DLL CPL_STDCALL
87 GDALRasterBandH hMaskBand,
88 OGRLayerH hOutLayer, int iPixValField,
89 char **papszOptions,
90 GDALProgressFunc pfnProgress,
91 void * pProgressArg );
92
93CPLErr CPL_DLL CPL_STDCALL
95 GDALRasterBandH hMaskBand,
96 OGRLayerH hOutLayer, int iPixValField,
97 char **papszOptions,
98 GDALProgressFunc pfnProgress,
99 void * pProgressArg );
100
101CPLErr CPL_DLL CPL_STDCALL
103 GDALRasterBandH hDstBand,
104 int nSizeThreshold, int nConnectedness,
105 char **papszOptions,
106 GDALProgressFunc pfnProgress,
107 void * pProgressArg );
108
109/*
110 * Warp Related.
111 */
112
113typedef int
114(*GDALTransformerFunc)( void *pTransformerArg,
115 int bDstToSrc, int nPointCount,
116 double *x, double *y, double *z, int *panSuccess );
117
119#define GDAL_GTI2_SIGNATURE "GTI2"
120
121typedef struct {
122 GByte abySignature[4];
123 const char *pszClassName;
124 GDALTransformerFunc pfnTransform;
125 void (*pfnCleanup)( void * pTransformerArg );
126 CPLXMLNode *(*pfnSerialize)( void * pTransformerArg );
127 void* (*pfnCreateSimilar)( void* pTransformerArg, double dfSrcRatioX, double dfSrcRatioY );
128} GDALTransformerInfo;
132void CPL_DLL GDALDestroyTransformer( void *pTransformerArg );
133int CPL_DLL GDALUseTransformer( void *pTransformerArg,
134 int bDstToSrc, int nPointCount,
135 double *x, double *y, double *z,
136 int *panSuccess );
137void* GDALCreateSimilarTransformer( void* psTransformerArg, double dfSrcRatioX, double dfSrcRatioY );
140/* High level transformer for going from image coordinates on one file
141 to image coordinates on another, potentially doing reprojection,
142 utilizing GCPs or using the geotransform. */
143
144void CPL_DLL *
145GDALCreateGenImgProjTransformer( GDALDatasetH hSrcDS, const char *pszSrcWKT,
146 GDALDatasetH hDstDS, const char *pszDstWKT,
147 int bGCPUseOK, double dfGCPErrorThreshold,
148 int nOrder );
149void CPL_DLL *
151 char **papszOptions );
152void CPL_DLL *
153GDALCreateGenImgProjTransformer3( const char *pszSrcWKT,
154 const double *padfSrcGeoTransform,
155 const char *pszDstWKT,
156 const double *padfDstGeoTransform );
157
158void CPL_DLL *
160 const double *padfSrcGeoTransform,
161 OGRSpatialReferenceH hDstSRS,
162 const double *padfDstGeoTransform,
163 const char* const *papszOptions );
164
166 const double * );
167void CPL_DLL GDALDestroyGenImgProjTransformer( void * );
168int CPL_DLL GDALGenImgProjTransform(
169 void *pTransformArg, int bDstToSrc, int nPointCount,
170 double *x, double *y, double *z, int *panSuccess );
171
172void GDALSetTransformerDstGeoTransform( void *, const double * );
173void GDALGetTransformerDstGeoTransform( void*, double* );
174
175/* Geo to geo reprojection transformer. */
176void CPL_DLL *
177GDALCreateReprojectionTransformer( const char *pszSrcWKT,
178 const char *pszDstWKT );
179void CPL_DLL *
181 OGRSpatialReferenceH hSrcSRS,
182 OGRSpatialReferenceH hDstSRS,
183 const char* const *papszOptions);
184void CPL_DLL GDALDestroyReprojectionTransformer( void * );
185int CPL_DLL GDALReprojectionTransform(
186 void *pTransformArg, int bDstToSrc, int nPointCount,
187 double *x, double *y, double *z, int *panSuccess );
188
189/* GCP based transformer ... forward is to georef coordinates */
190void CPL_DLL *
191GDALCreateGCPTransformer( int nGCPCount, const GDAL_GCP *pasGCPList,
192 int nReqOrder, int bReversed );
193
194/* GCP based transformer with refinement of the GCPs ... forward is to georef coordinates */
195void CPL_DLL *
196GDALCreateGCPRefineTransformer( int nGCPCount, const GDAL_GCP *pasGCPList,
197 int nReqOrder, int bReversed, double tolerance, int minimumGcps);
198
199void CPL_DLL GDALDestroyGCPTransformer( void *pTransformArg );
200int CPL_DLL GDALGCPTransform(
201 void *pTransformArg, int bDstToSrc, int nPointCount,
202 double *x, double *y, double *z, int *panSuccess );
203
204/* Thin Plate Spine transformer ... forward is to georef coordinates */
205
206void CPL_DLL *
207GDALCreateTPSTransformer( int nGCPCount, const GDAL_GCP *pasGCPList,
208 int bReversed );
209void CPL_DLL GDALDestroyTPSTransformer( void *pTransformArg );
210int CPL_DLL GDALTPSTransform(
211 void *pTransformArg, int bDstToSrc, int nPointCount,
212 double *x, double *y, double *z, int *panSuccess );
213
215char CPL_DLL ** RPCInfoToMD( GDALRPCInfo *psRPCInfo );
218/* RPC based transformer ... src is pixel/line/elev, dst is long/lat/elev */
219
220void CPL_DLL *
221GDALCreateRPCTransformer( GDALRPCInfo *psRPC, int bReversed,
222 double dfPixErrThreshold,
223 char **papszOptions );
224void CPL_DLL GDALDestroyRPCTransformer( void *pTransformArg );
225int CPL_DLL GDALRPCTransform(
226 void *pTransformArg, int bDstToSrc, int nPointCount,
227 double *x, double *y, double *z, int *panSuccess );
228
229/* Geolocation transformer */
230
231void CPL_DLL *
233 char **papszGeolocationInfo,
234 int bReversed );
235void CPL_DLL GDALDestroyGeoLocTransformer( void *pTransformArg );
236int CPL_DLL GDALGeoLocTransform(
237 void *pTransformArg, int bDstToSrc, int nPointCount,
238 double *x, double *y, double *z, int *panSuccess );
239
240/* Approximate transformer */
241void CPL_DLL *
243 void *pRawTransformerArg, double dfMaxError );
244void CPL_DLL GDALApproxTransformerOwnsSubtransformer( void *pCBData,
245 int bOwnFlag );
246void CPL_DLL GDALDestroyApproxTransformer( void *pApproxArg );
247int CPL_DLL GDALApproxTransform(
248 void *pTransformArg, int bDstToSrc, int nPointCount,
249 double *x, double *y, double *z, int *panSuccess );
250
251int CPL_DLL CPL_STDCALL
253 GDALDatasetH hDstDS,
254 int nBandCount, int *panBandList,
255 GDALTransformerFunc pfnTransform,
256 void *pTransformArg,
257 GDALProgressFunc pfnProgress,
258 void *pProgressArg,
259 char **papszWarpOptions );
260
261CPLErr CPL_DLL CPL_STDCALL
263 GDALTransformerFunc pfnTransformer,
264 void *pTransformArg,
265 double *padfGeoTransformOut,
266 int *pnPixels, int *pnLines );
267CPLErr CPL_DLL CPL_STDCALL
269 GDALTransformerFunc pfnTransformer,
270 void *pTransformArg,
271 double *padfGeoTransformOut,
272 int *pnPixels, int *pnLines,
273 double *padfExtents,
274 int nOptions );
275
277CPLXMLNode CPL_DLL *
278GDALSerializeTransformer( GDALTransformerFunc pfnFunc, void *pTransformArg );
279CPLErr CPL_DLL GDALDeserializeTransformer( CPLXMLNode *psTree,
280 GDALTransformerFunc *ppfnFunc,
281 void **ppTransformArg );
284CPLErr CPL_DLL
286 GDALRasterBandH hYBand,
287 GDALRasterBandH hZBand,
288 GDALTransformerFunc pfnTransformer,
289 void *pTransformArg,
290 GDALProgressFunc pfnProgress,
291 void *pProgressArg,
292 char **papszOptions );
293
294/* -------------------------------------------------------------------- */
295/* Contour Line Generation */
296/* -------------------------------------------------------------------- */
297
299typedef CPLErr (*GDALContourWriter)( double dfLevel, int nPoints,
300 double *padfX, double *padfY, void * );
301
304
306GDAL_CG_Create( int nWidth, int nHeight,
307 int bNoDataSet, double dfNoDataValue,
308 double dfContourInterval, double dfContourBase,
309 GDALContourWriter pfnWriter, void *pCBData );
311 double *padfScanline );
312void CPL_DLL GDAL_CG_Destroy( GDALContourGeneratorH hCG );
313
315typedef struct
316{
317 void *hLayer;
318
319 double adfGeoTransform[6];
320
321 int nElevField;
322 int nElevFieldMin;
323 int nElevFieldMax;
324 int nIDField;
325 int nNextID;
326} OGRContourWriterInfo;
327
328CPLErr CPL_DLL
329OGRContourWriter( double, int, double *, double *, void *pInfo );
332CPLErr CPL_DLL
334 double dfContourInterval, double dfContourBase,
335 int nFixedLevelCount, double *padfFixedLevels,
336 int bUseNoData, double dfNoDataValue,
337 void *hLayer, int iIDField, int iElevField,
338 GDALProgressFunc pfnProgress, void *pProgressArg );
339
340CPLErr CPL_DLL
341GDALContourGenerateEx( GDALRasterBandH hBand, void *hLayer,
342 CSLConstList options,
343 GDALProgressFunc pfnProgress, void *pProgressArg );
344
345/************************************************************************/
346/* Rasterizer API - geometries burned into GDAL raster. */
347/************************************************************************/
348
349CPLErr CPL_DLL
351 int nBandCount, int *panBandList,
352 int nGeomCount, OGRGeometryH *pahGeometries,
353 GDALTransformerFunc pfnTransformer,
354 void *pTransformArg,
355 double *padfGeomBurnValue,
356 char **papszOptions,
357 GDALProgressFunc pfnProgress,
358 void * pProgressArg );
359CPLErr CPL_DLL
361 int nBandCount, int *panBandList,
362 int nLayerCount, OGRLayerH *pahLayers,
363 GDALTransformerFunc pfnTransformer,
364 void *pTransformArg,
365 double *padfLayerBurnValues,
366 char **papszOptions,
367 GDALProgressFunc pfnProgress,
368 void *pProgressArg );
369
370CPLErr CPL_DLL
371GDALRasterizeLayersBuf( void *pData, int nBufXSize, int nBufYSize,
372 GDALDataType eBufType, int nPixelSpace, int nLineSpace,
373 int nLayerCount, OGRLayerH *pahLayers,
374 const char *pszDstProjection,
375 double *padfDstGeoTransform,
376 GDALTransformerFunc pfnTransformer,
377 void *pTransformArg, double dfBurnValue,
378 char **papszOptions, GDALProgressFunc pfnProgress,
379 void *pProgressArg );
380
381/************************************************************************/
382/* Gridding interface. */
383/************************************************************************/
384
402
404typedef struct
405{
407 double dfPower;
415 double dfRadius1;
417 double dfRadius2;
422 double dfAngle;
439
441typedef struct
442{
444 double dfPower;
446 double dfRadius;
449
466
468typedef struct
469{
471 double dfRadius1;
473 double dfRadius2;
478 double dfAngle;
488
490typedef struct
491{
493 double dfRadius1;
495 double dfRadius2;
500 double dfAngle;
504
506typedef struct
507{
509 double dfRadius1;
511 double dfRadius2;
516 double dfAngle;
526
528typedef struct
529{
535 double dfRadius;
539
540CPLErr CPL_DLL
542 const double *, const double *, const double *,
543 double, double, double, double,
544 GUInt32, GUInt32, GDALDataType, void *,
545 GDALProgressFunc, void *);
546
549
550GDALGridContext CPL_DLL*
551GDALGridContextCreate( GDALGridAlgorithm eAlgorithm, const void *poOptions,
552 GUInt32 nPoints,
553 const double *padfX, const double *padfY, const double *padfZ,
554 int bCallerWillKeepPointArraysAlive );
555
556void CPL_DLL GDALGridContextFree(GDALGridContext* psContext);
557
559 double dfXMin, double dfXMax, double dfYMin, double dfYMax,
560 GUInt32 nXSize, GUInt32 nYSize, GDALDataType eType, void *pData,
561 GDALProgressFunc pfnProgress, void *pProgressArg );
562
563GDAL_GCP CPL_DLL *
565 GDALDatasetH hSecondImage,
566 char **papszOptions,
567 int *pnGCPCount );
568
569/************************************************************************/
570/* Delaunay triangulation interface. */
571/************************************************************************/
572
574typedef struct
575{
576 int anVertexIdx[3];
577 int anNeighborIdx[3];
578 /* anNeighborIdx[k] is the triangle to the opposite side */
579 /* of the opposite segment of anVertexIdx[k] */
581
589typedef struct
590{
591 double dfMul1X;
592 double dfMul1Y;
593 double dfMul2X;
594 double dfMul2Y;
595 double dfCstX;
596 double dfCstY;
598
600typedef struct
601{
606
607int CPL_DLL GDALHasTriangulation(void);
608
610 const double* padfX,
611 const double* padfY);
613 GDALTriangulation* psDT,
614 const double* padfX,
615 const double* padfY);
617 const GDALTriangulation* psDT,
618 int nFacetIdx,
619 double dfX,
620 double dfY,
621 double* pdfL1,
622 double* pdfL2,
623 double* pdfL3);
625 double dfX,
626 double dfY,
627 int* panOutputFacetIdx );
629 int nFacetIdx,
630 double dfX,
631 double dfY,
632 int* panOutputFacetIdx );
633void CPL_DLL GDALTriangulationFree(GDALTriangulation* psDT);
634
636// GDAL internal use only
637void GDALTriangulationTerminate(void);
641 const char* pszProj4Geoidgrids,
642 int* pbError );
643
645 GDALDatasetH hGridDataset,
646 int bInverse,
647 double dfSrcUnitToMeter,
648 double dfDstUnitToMeter,
649 const char* const* papszOptions );
650
652
653#endif /* ndef GDAL_ALG_H_INCLUDED */
CPLErr
Error category.
Definition: cpl_error.h:53
Definitions for CPL mini XML Parser/Serializer.
#define CPL_C_END
Macro to end a block of C symbols.
Definition: cpl_port.h:339
#define CPL_C_START
Macro to start a block of C symbols.
Definition: cpl_port.h:337
unsigned int GUInt32
Unsigned int32 type.
Definition: cpl_port.h:207
char ** CSLConstList
Type of a constant null-terminated list of nul terminated strings.
Definition: cpl_port.h:1194
unsigned char GByte
Unsigned byte type.
Definition: cpl_port.h:215
Public (C callable) GDAL entry points.
GDALDataType
Definition: gdal.h:60
void * GDALDatasetH
Opaque type used for the C bindings of the C++ GDALDataset class.
Definition: gdal.h:255
void * GDALRasterBandH
Opaque type used for the C bindings of the C++ GDALRasterBand class.
Definition: gdal.h:258
void * GDALColorTableH
Opaque type used for the C bindings of the C++ GDALColorTable class.
Definition: gdal.h:264
CPLErr GDALContourGenerateEx(GDALRasterBandH hBand, void *hLayer, CSLConstList options, GDALProgressFunc pfnProgress, void *pProgressArg)
Create vector contours from raster DEM.
Definition: contour.cpp:523
void * GDALCreateGeoLocTransformer(GDALDatasetH hBaseDS, char **papszGeolocationInfo, int bReversed)
Create GeoLocation transformer.
Definition: gdalgeoloc.cpp:650
struct GDALGridContext GDALGridContext
Grid context opaque type.
Definition: gdal_alg.h:548
CPLErr GDALFillNodata(GDALRasterBandH hTargetBand, GDALRasterBandH hMaskBand, double dfMaxSearchDist, int bDeprecatedOption, int nSmoothingIterations, char **papszOptions, GDALProgressFunc pfnProgress, void *pProgressArg)
Fill selected raster regions by interpolation from the edges.
Definition: rasterfill.cpp:412
int GDALGenImgProjTransform(void *pTransformArg, int bDstToSrc, int nPointCount, double *x, double *y, double *z, int *panSuccess)
Perform general image reprojection transformation.
Definition: gdaltransformer.cpp:2297
GDALTriangulation * GDALTriangulationCreateDelaunay(int nPoints, const double *padfX, const double *padfY)
Computes a Delaunay triangulation of the passed points.
Definition: delaunay.c:118
GDALContourGeneratorH GDAL_CG_Create(int nWidth, int nHeight, int bNoDataSet, double dfNoDataValue, double dfContourInterval, double dfContourBase, GDALContourWriter pfnWriter, void *pCBData)
Create contour generator.
Definition: contour.cpp:711
void * GDALCreateGCPRefineTransformer(int nGCPCount, const GDAL_GCP *pasGCPList, int nReqOrder, int bReversed, double tolerance, int minimumGcps)
Create GCP based polynomial transformer, with a tolerance threshold to discard GCPs that transform ba...
Definition: gdal_crs.cpp:341
void * GDALCreateReprojectionTransformer(const char *pszSrcWKT, const char *pszDstWKT)
Create reprojection transformer.
Definition: gdaltransformer.cpp:2722
void * GDALCreateReprojectionTransformerEx(OGRSpatialReferenceH hSrcSRS, OGRSpatialReferenceH hDstSRS, const char *const *papszOptions)
Create reprojection transformer.
Definition: gdaltransformer.cpp:2789
void GDALSetGenImgProjTransformerDstGeoTransform(void *, const double *)
Set GenImgProj output geotransform.
Definition: gdaltransformer.cpp:2228
void GDALDestroyGeoLocTransformer(void *pTransformArg)
Destroy GeoLocation transformer.
Definition: gdalgeoloc.cpp:823
void GDALGridContextFree(GDALGridContext *psContext)
Free a context used created by GDALGridContextCreate()
Definition: gdalgrid.cpp:2129
void GDALTriangulationFree(GDALTriangulation *psDT)
Free a triangulation.
Definition: delaunay.c:275
void * GDALCreateGenImgProjTransformer(GDALDatasetH hSrcDS, const char *pszSrcWKT, GDALDatasetH hDstDS, const char *pszDstWKT, int bGCPUseOK, double dfGCPErrorThreshold, int nOrder)
Create image to image transformer.
Definition: gdaltransformer.cpp:1070
CPLErr GDALGridCreate(GDALGridAlgorithm, const void *, GUInt32, const double *, const double *, const double *, double, double, double, double, GUInt32, GUInt32, GDALDataType, void *, GDALProgressFunc, void *)
Create regular grid from the scattered data.
Definition: gdalgrid.cpp:2422
int GDALDitherRGB2PCT(GDALRasterBandH hRed, GDALRasterBandH hGreen, GDALRasterBandH hBlue, GDALRasterBandH hTarget, GDALColorTableH hColorTable, GDALProgressFunc pfnProgress, void *pProgressArg)
24bit to 8bit conversion with dithering.
Definition: gdaldither.cpp:141
int GDALTriangulationComputeBarycentricCoordinates(const GDALTriangulation *psDT, int nFacetIdx, double dfX, double dfY, double *pdfL1, double *pdfL2, double *pdfL3)
Computes the barycentric coordinates of a point.
Definition: delaunay.c:378
CPLErr GDALSieveFilter(GDALRasterBandH hSrcBand, GDALRasterBandH hMaskBand, GDALRasterBandH hDstBand, int nSizeThreshold, int nConnectedness, char **papszOptions, GDALProgressFunc pfnProgress, void *pProgressArg)
Removes small raster polygons.
Definition: gdalsievefilter.cpp:200
GDALDatasetH GDALApplyVerticalShiftGrid(GDALDatasetH hSrcDataset, GDALDatasetH hGridDataset, int bInverse, double dfSrcUnitToMeter, double dfDstUnitToMeter, const char *const *papszOptions)
Apply a vertical shift grid to a source (DEM typically) dataset.
Definition: gdalapplyverticalshiftgrid.cpp:374
CPLErr GDALPolygonize(GDALRasterBandH hSrcBand, GDALRasterBandH hMaskBand, OGRLayerH hOutLayer, int iPixValField, char **papszOptions, GDALProgressFunc pfnProgress, void *pProgressArg)
Create polygon coverage from raster data.
Definition: polygonize.cpp:826
void GDALDestroyGCPTransformer(void *pTransformArg)
Destroy GCP transformer.
Definition: gdal_crs.cpp:368
CPLErr GDALRasterizeGeometries(GDALDatasetH hDS, int nBandCount, int *panBandList, int nGeomCount, OGRGeometryH *pahGeometries, GDALTransformerFunc pfnTransformer, void *pTransformArg, double *padfGeomBurnValue, char **papszOptions, GDALProgressFunc pfnProgress, void *pProgressArg)
Burn geometries into raster.
Definition: gdalrasterize.cpp:722
GDALGridAlgorithm
Gridding Algorithms.
Definition: gdal_alg.h:386
@ GGA_MetricMinimum
Definition: gdal_alg.h:390
@ GGA_InverseDistanceToAPowerNearestNeighbor
Definition: gdal_alg.h:400
@ GGA_InverseDistanceToAPower
Definition: gdal_alg.h:387
@ GGA_MetricAverageDistancePts
Definition: gdal_alg.h:396
@ GGA_MetricMaximum
Definition: gdal_alg.h:391
@ GGA_NearestNeighbor
Definition: gdal_alg.h:389
@ GGA_MetricAverageDistance
Definition: gdal_alg.h:394
@ GGA_MovingAverage
Definition: gdal_alg.h:388
@ GGA_MetricCount
Definition: gdal_alg.h:393
@ GGA_MetricRange
Definition: gdal_alg.h:392
@ GGA_Linear
Definition: gdal_alg.h:398
void GDALDestroyGenImgProjTransformer(void *)
GenImgProjTransformer deallocator.
Definition: gdaltransformer.cpp:2259
int GDALRPCTransform(void *pTransformArg, int bDstToSrc, int nPointCount, double *x, double *y, double *z, int *panSuccess)
RPC transform.
Definition: gdal_rpc.cpp:2031
int GDALComputeMedianCutPCT(GDALRasterBandH hRed, GDALRasterBandH hGreen, GDALRasterBandH hBlue, int(*pfnIncludePixel)(int, int, void *), int nColors, GDALColorTableH hColorTable, GDALProgressFunc pfnProgress, void *pProgressArg)
Compute optimal PCT for RGB image.
Definition: gdalmediancut.cpp:147
int GDALTriangulationFindFacetBruteForce(const GDALTriangulation *psDT, double dfX, double dfY, int *panOutputFacetIdx)
Returns the index of the triangle that contains the point by iterating over all triangles.
Definition: delaunay.c:426
GDALDatasetH GDALOpenVerticalShiftGrid(const char *pszProj4Geoidgrids, int *pbError)
Load proj.4 geoidgrids as GDAL dataset.
Definition: gdalapplyverticalshiftgrid.cpp:615
void GDAL_CG_Destroy(GDALContourGeneratorH hCG)
Destroy contour generator.
Definition: contour.cpp:745
void GDALApproxTransformerOwnsSubtransformer(void *pCBData, int bOwnFlag)
Set bOwnSubtransformer flag.
Definition: gdaltransformer.cpp:3273
int GDALApproxTransform(void *pTransformArg, int bDstToSrc, int nPointCount, double *x, double *y, double *z, int *panSuccess)
Perform approximate transformation.
Definition: gdaltransformer.cpp:3612
void * GDALCreateGenImgProjTransformer4(OGRSpatialReferenceH hSrcSRS, const double *padfSrcGeoTransform, OGRSpatialReferenceH hDstSRS, const double *padfDstGeoTransform, const char *const *papszOptions)
Create image to image transformer.
Definition: gdaltransformer.cpp:2119
CPLErr GDALFPolygonize(GDALRasterBandH hSrcBand, GDALRasterBandH hMaskBand, OGRLayerH hOutLayer, int iPixValField, char **papszOptions, GDALProgressFunc pfnProgress, void *pProgressArg)
Create polygon coverage from raster data.
Definition: polygonize.cpp:904
CPLErr GDALGridContextProcess(GDALGridContext *psContext, double dfXMin, double dfXMax, double dfYMin, double dfYMax, GUInt32 nXSize, GUInt32 nYSize, GDALDataType eType, void *pData, GDALProgressFunc pfnProgress, void *pProgressArg)
Do the gridding of a window of a raster.
Definition: gdalgrid.cpp:2183
int(* GDALTransformerFunc)(void *pTransformerArg, int bDstToSrc, int nPointCount, double *x, double *y, double *z, int *panSuccess)
Definition: gdal_alg.h:114
GDALGridContext * GDALGridContextCreate(GDALGridAlgorithm eAlgorithm, const void *poOptions, GUInt32 nPoints, const double *padfX, const double *padfY, const double *padfZ, int bCallerWillKeepPointArraysAlive)
Creates a context to do regular gridding from the scattered data.
Definition: gdalgrid.cpp:1733
CPLErr GDALSuggestedWarpOutput(GDALDatasetH hSrcDS, GDALTransformerFunc pfnTransformer, void *pTransformArg, double *padfGeoTransformOut, int *pnPixels, int *pnLines)
Suggest output file size.
Definition: gdaltransformer.cpp:177
CPLErr GDALComputeProximity(GDALRasterBandH hSrcBand, GDALRasterBandH hProximityBand, char **papszOptions, GDALProgressFunc pfnProgress, void *pProgressArg)
Compute the proximity of all pixels in the image to a set of pixels in the source image.
Definition: gdalproximity.cpp:112
CPLErr GDALTransformGeolocations(GDALRasterBandH hXBand, GDALRasterBandH hYBand, GDALRasterBandH hZBand, GDALTransformerFunc pfnTransformer, void *pTransformArg, GDALProgressFunc pfnProgress, void *pProgressArg, char **papszOptions)
Transform locations held in bands.
Definition: gdaltransformgeolocs.cpp:68
void GDALGetTransformerDstGeoTransform(void *, double *)
Get ApproxTransformer or GenImgProj output geotransform.
Definition: gdaltransformer.cpp:4259
CPLErr GDALSuggestedWarpOutput2(GDALDatasetH hSrcDS, GDALTransformerFunc pfnTransformer, void *pTransformArg, double *padfGeoTransformOut, int *pnPixels, int *pnLines, double *padfExtents, int nOptions)
Suggest output file size.
Definition: gdaltransformer.cpp:354
GDAL_GCP * GDALComputeMatchingPoints(GDALDatasetH hFirstImage, GDALDatasetH hSecondImage, char **papszOptions, int *pnGCPCount)
GDALComputeMatchingPoints.
Definition: gdalmatching.cpp:188
CPLErr GDAL_CG_FeedLine(GDALContourGeneratorH hCG, double *padfScanline)
Feed a line to the contour generator.
Definition: contour.cpp:733
int GDALTriangulationComputeBarycentricCoefficients(GDALTriangulation *psDT, const double *padfX, const double *padfY)
Computes barycentric coefficients for each triangles of the triangulation.
Definition: delaunay.c:301
int GDALGCPTransform(void *pTransformArg, int bDstToSrc, int nPointCount, double *x, double *y, double *z, int *panSuccess)
Transforms point based on GCP derived polynomial model.
Definition: gdal_crs.cpp:410
void * GDALCreateRPCTransformer(GDALRPCInfo *psRPC, int bReversed, double dfPixErrThreshold, char **papszOptions)
Create an RPC based transformer.
Definition: gdal_rpc.cpp:776
void GDALDestroyApproxTransformer(void *pApproxArg)
Cleanup approximate transformer.
Definition: gdaltransformer.cpp:3294
void * GDALCreateGenImgProjTransformer3(const char *pszSrcWKT, const double *padfSrcGeoTransform, const char *pszDstWKT, const double *padfDstGeoTransform)
Create image to image transformer.
Definition: gdaltransformer.cpp:2064
int GDALHasTriangulation(void)
Returns if GDAL is built with Delaunay triangulation support.
Definition: delaunay.c:95
CPLErr(* GDALContourWriter)(double dfLevel, int nPoints, double *padfX, double *padfY, void *)
Contour writer callback type.
Definition: gdal_alg.h:299
CPLErr GDALRasterizeLayersBuf(void *pData, int nBufXSize, int nBufYSize, GDALDataType eBufType, int nPixelSpace, int nLineSpace, int nLayerCount, OGRLayerH *pahLayers, const char *pszDstProjection, double *padfDstGeoTransform, GDALTransformerFunc pfnTransformer, void *pTransformArg, double dfBurnValue, char **papszOptions, GDALProgressFunc pfnProgress, void *pProgressArg)
Burn geometries from the specified list of layer into raster.
Definition: gdalrasterize.cpp:1540
int GDALTPSTransform(void *pTransformArg, int bDstToSrc, int nPointCount, double *x, double *y, double *z, int *panSuccess)
Transforms point based on GCP derived polynomial model.
Definition: gdal_tps.cpp:358
int GDALChecksumImage(GDALRasterBandH hBand, int nXOff, int nYOff, int nXSize, int nYSize)
Compute checksum for image region.
Definition: gdalchecksum.cpp:66
void GDALSetTransformerDstGeoTransform(void *, const double *)
Set ApproxTransformer or GenImgProj output geotransform.
Definition: gdaltransformer.cpp:4235
void * GDALContourGeneratorH
Contour generator opaque type.
Definition: gdal_alg.h:303
CPLErr GDALContourGenerate(GDALRasterBandH hBand, double dfContourInterval, double dfContourBase, int nFixedLevelCount, double *padfFixedLevels, int bUseNoData, double dfNoDataValue, void *hLayer, int iIDField, int iElevField, GDALProgressFunc pfnProgress, void *pProgressArg)
Create vector contours from raster DEM.
Definition: contour.cpp:310
void GDALDestroyRPCTransformer(void *pTransformArg)
Destroy RPC tranformer.
Definition: gdal_rpc.cpp:1054
CPLErr GDALRasterizeLayers(GDALDatasetH hDS, int nBandCount, int *panBandList, int nLayerCount, OGRLayerH *pahLayers, GDALTransformerFunc pfnTransformer, void *pTransformArg, double *padfLayerBurnValues, char **papszOptions, GDALProgressFunc pfnProgress, void *pProgressArg)
Burn geometries from the specified list of layers into raster.
Definition: gdalrasterize.cpp:1147
void * GDALCreateGenImgProjTransformer2(GDALDatasetH hSrcDS, GDALDatasetH hDstDS, char **papszOptions)
Create image to image transformer.
Definition: gdaltransformer.cpp:1470
int GDALSimpleImageWarp(GDALDatasetH hSrcDS, GDALDatasetH hDstDS, int nBandCount, int *panBandList, GDALTransformerFunc pfnTransform, void *pTransformArg, GDALProgressFunc pfnProgress, void *pProgressArg, char **papszWarpOptions)
Perform simple image warp.
Definition: gdalsimplewarp.cpp:230
void GDALDestroyReprojectionTransformer(void *)
Destroy reprojection transformation.
Definition: gdaltransformer.cpp:2891
void * GDALCreateGCPTransformer(int nGCPCount, const GDAL_GCP *pasGCPList, int nReqOrder, int bReversed)
Create GCP based polynomial transformer.
Definition: gdal_crs.cpp:331
int GDALReprojectionTransform(void *pTransformArg, int bDstToSrc, int nPointCount, double *x, double *y, double *z, int *panSuccess)
Perform reprojection transformation.
Definition: gdaltransformer.cpp:2924
void * GDALCreateApproxTransformer(GDALTransformerFunc pfnRawTransformer, void *pRawTransformerArg, double dfMaxError)
Create an approximating transformer.
Definition: gdaltransformer.cpp:3230
int GDALGeoLocTransform(void *pTransformArg, int bDstToSrc, int nPointCount, double *x, double *y, double *z, int *panSuccess)
Use GeoLocation transformer.
Definition: gdalgeoloc.cpp:854
int GDALTriangulationFindFacetDirected(const GDALTriangulation *psDT, int nFacetIdx, double dfX, double dfY, int *panOutputFacetIdx)
Returns the index of the triangle that contains the point by walking in the triangulation.
Definition: delaunay.c:521
void * GDALCreateTPSTransformer(int nGCPCount, const GDAL_GCP *pasGCPList, int bReversed)
Create Thin Plate Spline transformer from GCPs.
Definition: gdal_tps.cpp:144
void GDALDestroyTPSTransformer(void *pTransformArg)
Destroy TPS transformer.
Definition: gdal_tps.cpp:313
C API and defines for OGRFeature, OGRGeometry, and OGRDataSource related classes.
void * OGRGeometryH
Opaque type for a geometry.
Definition: ogr_api.h:60
void * OGRSpatialReferenceH
Opaque type for a spatial reference system.
Definition: ogr_api.h:74
void * OGRLayerH
Opaque type for a layer (OGRLayer)
Definition: ogr_api.h:509
Document node structure.
Definition: cpl_minixml.h:67
Data metrics method control options.
Definition: gdal_alg.h:507
double dfRadius2
Definition: gdal_alg.h:511
double dfAngle
Definition: gdal_alg.h:516
double dfNoDataValue
Definition: gdal_alg.h:524
double dfRadius1
Definition: gdal_alg.h:509
GUInt32 nMinPoints
Definition: gdal_alg.h:522
Inverse distance to a power, with nearest neighbour search, control options.
Definition: gdal_alg.h:442
GUInt32 nMinPoints
Definition: gdal_alg.h:462
double dfNoDataValue
Definition: gdal_alg.h:464
GUInt32 nMaxPoints
Definition: gdal_alg.h:456
double dfSmoothing
Definition: gdal_alg.h:448
Inverse distance to a power method control options.
Definition: gdal_alg.h:405
GUInt32 nMinPoints
Definition: gdal_alg.h:435
double dfPower
Definition: gdal_alg.h:407
double dfRadius2
Definition: gdal_alg.h:417
double dfNoDataValue
Definition: gdal_alg.h:437
double dfSmoothing
Definition: gdal_alg.h:409
double dfAngle
Definition: gdal_alg.h:422
double dfRadius1
Definition: gdal_alg.h:415
GUInt32 nMaxPoints
Definition: gdal_alg.h:429
double dfAnisotropyRatio
Definition: gdal_alg.h:411
double dfAnisotropyAngle
Definition: gdal_alg.h:413
Linear method control options.
Definition: gdal_alg.h:529
double dfNoDataValue
Definition: gdal_alg.h:537
double dfRadius
Definition: gdal_alg.h:535
Moving average method control options.
Definition: gdal_alg.h:469
double dfRadius1
Definition: gdal_alg.h:471
double dfNoDataValue
Definition: gdal_alg.h:486
double dfAngle
Definition: gdal_alg.h:478
GUInt32 nMinPoints
Definition: gdal_alg.h:484
double dfRadius2
Definition: gdal_alg.h:473
Nearest neighbor method control options.
Definition: gdal_alg.h:491
double dfNoDataValue
Definition: gdal_alg.h:502
double dfRadius2
Definition: gdal_alg.h:495
double dfAngle
Definition: gdal_alg.h:500
double dfRadius1
Definition: gdal_alg.h:493
Structure to store Rational Polynomial Coefficients / Rigorous Projection Model.
Definition: gdal.h:1028
Triangle barycentric coefficients.
Definition: gdal_alg.h:590
double dfCstY
dfCstY
Definition: gdal_alg.h:596
double dfCstX
dfCstX
Definition: gdal_alg.h:595
double dfMul2Y
dfMul2Y
Definition: gdal_alg.h:594
double dfMul1Y
dfMul1Y
Definition: gdal_alg.h:592
double dfMul2X
dfMul2X
Definition: gdal_alg.h:593
double dfMul1X
dfMul1X
Definition: gdal_alg.h:591
Triangle fact.
Definition: gdal_alg.h:575
Triangulation structure.
Definition: gdal_alg.h:601
GDALTriBarycentricCoefficients * pasFacetCoefficients
arra of nFacets barycentric coefficients
Definition: gdal_alg.h:604
int nFacets
number of facets
Definition: gdal_alg.h:602
GDALTriFacet * pasFacets
array of nFacets facets
Definition: gdal_alg.h:603
Ground Control Point.
Definition: gdal.h:564

Generated for GDAL by doxygen 1.9.4.