FortranGIS Version 3.0
gdal Module Reference

Fortran 2003 interface to the gdal http://www.gdal.org/ library. More...

Data Types

interface  gdalapplygeotransform_f
 Interface to a Fortran version of gdalapplygeotransform working on scalars, 1-d, 2-d and 3-d arrays. More...
 
interface  gdaldatasetrasterio_f
 Simplified Fortran generic interface to the gdaldatasetrasterio C function. More...
 
interface  gdalmajorobjecth_new
 Fortran interface for formally converting a dataset, rasterband or driver opaque object into a generic gdal object of type gdalmajorobjecth to be used in some methods such as GDALGetMetadata. More...
 
interface  gdalrasterio_f
 Simplified Fortran generic interface to the gdalrasterio C function. More...
 

Functions/Subroutines

subroutine gdaldatasetbbsize_f (hds, bbxmin, bbymin, bbxmax, bbymax, nx, ny, offsetx, offsety, xmin, ymin, xmax, ymax)
 Determine the size of a dataset subset.
 
subroutine gdaldatasetsimpleread_f (hds, bbxmin, bbymin, bbxmax, bbymax, pbuffer, xmin, ymin, xmax, ymax)
 Even more simplified method for importing data from a dataset within a bounding box specified in georeferenced coordinates.
 
subroutine gdalrastersimpleread_f (hband, bbxmin, bbymin, bbxmax, bbymax, pbuffer, xmin, ymin, xmax, ymax)
 Even more simplified method for importing data from a raster band within a bounding box specified in georeferenced coordinates.
 

Variables

integer(kind=c_int), parameter gdt_unknown = 0
 constant defining the native data type of a dataset data: unknown
 
integer(kind=c_int), parameter gdt_byte = 1
 byte, in Fortran it can be declared as INTEGER(kind=C_INT_8_T)
 
integer(kind=c_int), parameter gdt_uint16 = 2
 unsigned 16 bit integer, it should be avoided in Fortran and translated into a signed type
 
integer(kind=c_int), parameter gdt_int16 = 3
 signed 16 bit integer, in Fortran it can be declared as INTEGER(kind=C_INT_16_T)
 
integer(kind=c_int), parameter gdt_uint32 = 4
 unsigned 32 bit integer, it should be avoided in Fortran and translated into a signed type
 
integer(kind=c_int), parameter gdt_int32 = 5
 signed 32 bit integer, in Fortran it can be declared as INTEGER(kind=C_INT)
 
integer(kind=c_int), parameter gdt_float32 = 6
 32 bit floating point real, in Fortran it can be declared as REAL(kind=C_FLOAT)
 
integer(kind=c_int), parameter gdt_float64 = 7
 64 bit floating point real, in Fortran it can be declared as REAL(kind=C_DOUBLE)
 
integer(kind=c_int), parameter gdt_cint16 = 8
 16 bit integer complex, it should be avoided in Fortran and translated into a floating point type
 
integer(kind=c_int), parameter gdt_cint32 = 9
 32 bit integer complex, it should be avoided in Fortran and translated into a floating point type
 
integer(kind=c_int), parameter gdt_cfloat32 = 10
 32 bit (*2) floating point complex, in Fortran it can be declared as COMPLEX(kind=C_FLOAT_COMPLEX)
 
integer(kind=c_int), parameter gdt_cfloat64 = 11
 64 bit (*2) floating point complex, in Fortran it can be declared as COMPLEX(kind=C_DOUBLE_COMPLEX)
 
integer(kind=c_int), parameter ga_readonly = 0
 access type for opening a file: read only
 
integer(kind=c_int), parameter ga_update = 1
 update
 
integer(kind=c_int), parameter gf_read = 0
 operation to be performed on a dataset: read
 
integer(kind=c_int), parameter gf_write = 1
 write
 

Detailed Description

Fortran 2003 interface to the gdal http://www.gdal.org/ library.

This module includes the interface to most of the public gdal C API plus a more "Fortranic" version of some functions.

The following functions are directly interfaced to their corresponding C version, so they are undocumented here, please refer to the original gdal C API documentation, e.g. at the address http://www.gdal.org/gdal_8h.html , for their use:

  • GDALDataTypeUnion() -> FUNCTION gdaldatatypeunion()
  • GDALGetDataTypeSize() -> FUNCTION gdalgetdatatypesize()
  • GDALDataTypeIsComplex() -> FUNCTION gdaldatatypeiscomplex()
  • GDALGetDataTypeName() -> FUNCTION gdalgetdatatypename()
  • GDALGetDataTypeByName() -> FUNCTION gdalgetdatatypebyname()
  • GDALGetPaletteInterpretationName() -> FUNCTION gdalgetpaletteinterpretationname()
  • GDALGetColorInterpretationName() -> FUNCTION gdalgetcolorinterpretationname()
  • GDALGetColorInterpretationByName() -> FUNCTION gdalgetcolorinterpretationbyname()
  • GDALAllRegister() -> SUBROUTINE gdalallregister()
  • GDALCreate() -> FUNCTION gdalcreate()
  • GDALCreateCopy() -> FUNCTION gdalcreatecopy()
  • GDALIdentifyDriver() -> FUNCTION gdalidentifydriver()
  • GDALOpen() -> FUNCTION gdalopen()
  • GDALOpenShared() -> FUNCTION gdalopenshared()
  • GDALGetDriverByName() -> FUNCTION gdalgetdriverbyname()
  • GDALGetDriverCount() -> FUNCTION gdalgetdrivercount()
  • GDALGetDriver() -> FUNCTION gdalgetdriver()
  • GDALDestroyDriver() -> SUBROUTINE gdaldestroydriver()
  • GDALRegisterDriver() -> FUNCTION gdalregisterdriver()
  • GDALDeregisterDriver() -> SUBROUTINE gdalderegisterdriver()
  • GDALDestroyDriverManager() -> SUBROUTINE gdaldestroydrivermanager()
  • GDALDeleteDataset() -> FUNCTION gdaldeletedataset()
  • GDALRenameDataset() -> FUNCTION gdalrenamedataset()
  • GDALCopyDatasetFiles() -> FUNCTION gdalcopydatasetfiles()
  • GDALValidateCreationOptions() -> FUNCTION gdalvalidatecreationoptions()
  • GDALGetDriverShortName() -> FUNCTION gdalgetdrivershortname()
  • GDALGetDriverLongName() -> FUNCTION gdalgetdriverlongname()
  • GDALGetDriverHelpTopic() -> FUNCTION gdalgetdriverhelptopic()
  • GDALGetDriverCreationOptionList() -> FUNCTION gdalgetdrivercreationoptionlist()
  • GDALInvGeoTransform() -> FUNCTION gdalinvgeotransform()
  • GDALGetMetadata() -> FUNCTION gdalgetmetadata()
  • GDALSetMetadata() -> FUNCTION gdalsetmetadata()
  • GDALGetMetadataItem() -> FUNCTION gdalgetmetadataitem()
  • GDALSetMetadataItem() -> FUNCTION gdalsetmetadataitem()
  • GDALGetDescription() -> FUNCTION gdalgetdescription()
  • GDALSetDescription() -> SUBROUTINE gdalsetdescription()
  • GDALGetDatasetDriver() -> FUNCTION gdalgetdatasetdriver()
  • GDALGetFileList() -> FUNCTION gdalgetfilelist()
  • GDALClose() -> SUBROUTINE gdalclose()
  • GDALGetRasterXSize() -> FUNCTION gdalgetrasterxsize()
  • GDALGetRasterYSize() -> FUNCTION gdalgetrasterysize()
  • GDALGetRasterCount() -> FUNCTION gdalgetrastercount()
  • GDALGetRasterBand() -> FUNCTION gdalgetrasterband()
  • GDALAddBand() -> FUNCTION gdaladdband()
  • GDALBeginAsyncReader() -> FUNCTION gdalbeginasyncreader()
  • GDALEndAsyncReader() -> SUBROUTINE gdalendasyncreader()
  • GDALDatasetRasterIO() -> FUNCTION gdaldatasetrasterio()
  • GDALDatasetAdviseRead() -> FUNCTION gdaldatasetadviseread()
  • GDALGetProjectionRef() -> FUNCTION gdalgetprojectionref()
  • GDALSetProjection() -> FUNCTION gdalsetprojection()
  • GDALGetGeoTransform() -> FUNCTION gdalgetgeotransform()
  • GDALSetGeoTransform() -> FUNCTION gdalsetgeotransform()
  • GDALGetGCPCount() -> FUNCTION gdalgetgcpcount()
  • GDALGetGCPProjection() -> FUNCTION gdalgetgcpprojection()
  • GDALGetInternalHandle() -> FUNCTION gdalgetinternalhandle()
  • GDALReferenceDataset() -> FUNCTION gdalreferencedataset()
  • GDALDereferenceDataset() -> FUNCTION gdaldereferencedataset()
  • GDALBuildOverviews() -> FUNCTION gdalbuildoverviews()
  • GDALGetOpenDatasets() -> SUBROUTINE gdalgetopendatasets()
  • GDALGetAccess() -> FUNCTION gdalgetaccess()
  • GDALFlushCache() -> SUBROUTINE gdalflushcache()
  • GDALCreateDatasetMaskBand() -> FUNCTION gdalcreatedatasetmaskband()
  • GDALDatasetCopyWholeRaster() -> FUNCTION gdaldatasetcopywholeraster()
  • GDALRasterBandCopyWholeRaster() -> FUNCTION gdalrasterbandcopywholeraster()
  • GDALRegenerateOverviews() -> FUNCTION gdalregenerateoverviews()
  • GDALGetRasterDataType() -> FUNCTION gdalgetrasterdatatype()
  • GDALRasterAdviseRead() -> FUNCTION gdalrasteradviseread()
  • GDALRasterIO() -> FUNCTION gdalrasterio()
  • GDALGetRasterBandXSize() -> FUNCTION gdalgetrasterbandxsize()
  • GDALGetRasterBandYSize() -> FUNCTION gdalgetrasterbandysize()
  • GDALGetRasterAccess() -> FUNCTION gdalgetrasteraccess()
  • GDALGetBandNumber() -> FUNCTION gdalgetbandnumber()
  • GDALGetBandDataset() -> FUNCTION gdalgetbanddataset()
  • GDALGetRasterColorInterpretation() -> FUNCTION gdalgetrastercolorinterpretation()
  • GDALSetRasterColorInterpretation() -> FUNCTION gdalsetrastercolorinterpretation()
  • GDALGetRasterColorTable() -> FUNCTION gdalgetrastercolortable()
  • GDALSetRasterColorTable() -> FUNCTION gdalsetrastercolortable()
  • GDALHasArbitraryOverviews() -> FUNCTION gdalhasarbitraryoverviews()
  • GDALGetOverviewCount() -> FUNCTION gdalgetoverviewcount()
  • GDALGetOverview() -> FUNCTION gdalgetoverview()
  • GDALGetRasterNoDataValue() -> FUNCTION gdalgetrasternodatavalue()
  • GDALSetRasterNoDataValue() -> FUNCTION gdalsetrasternodatavalue()
  • GDALGetRasterCategoryNames() -> FUNCTION gdalgetrastercategorynames()
  • GDALGetRasterMinimum() -> FUNCTION gdalgetrasterminimum()
  • GDALGetRasterMaximum() -> FUNCTION gdalgetrastermaximum()
  • GDALGetRasterStatistics() -> FUNCTION gdalgetrasterstatistics()
  • GDALComputeRasterStatistics() -> FUNCTION gdalcomputerasterstatistics()
  • GDALSetRasterStatistics() -> FUNCTION gdalsetrasterstatistics()
  • GDALGetRasterUnitType() -> FUNCTION gdalgetrasterunittype()
  • GDALSetRasterUnitType() -> FUNCTION gdalsetrasterunittype()
  • GDALGetRasterOffset() -> FUNCTION gdalgetrasteroffset()
  • GDALSetRasterOffset() -> FUNCTION gdalsetrasteroffset()
  • GDALGetRasterScale() -> FUNCTION gdalgetrasterscale()
  • GDALSetRasterScale() -> FUNCTION gdalsetrasterscale()
  • GDALComputeRasterMinMax() -> SUBROUTINE gdalcomputerasterminmax()
  • GDALFlushRasterCache() -> FUNCTION gdalflushrastercache()
  • GDALGetRasterHistogram() -> FUNCTION gdalgetrasterhistogram()
  • GDALGetDefaultHistogram() -> FUNCTION gdalgetdefaulthistogram()
  • GDALSetDefaultHistogram() -> FUNCTION gdalsetdefaulthistogram()
  • GDALGetRandomRasterSample() -> FUNCTION gdalgetrandomrastersample()
  • GDALGetRasterSampleOverview() -> FUNCTION gdalgetrastersampleoverview()
  • GDALFillRaster() -> FUNCTION gdalfillraster()
  • GDALComputeBandStats() -> FUNCTION gdalcomputebandstats()
  • GDALOverviewMagnitudeCorrection() -> FUNCTION gdaloverviewmagnitudecorrection()
  • GDALGetDefaultRAT() -> FUNCTION gdalgetdefaultrat()
  • GDALSetDefaultRAT() -> FUNCTION gdalsetdefaultrat()
  • GDALAddDerivedBandPixelFunc() -> FUNCTION gdaladdderivedbandpixelfunc()
  • GDALGetMaskBand() -> FUNCTION gdalgetmaskband()
  • GDALGetMaskFlags() -> FUNCTION gdalgetmaskflags()
  • GDALCreateMaskBand() -> FUNCTION gdalcreatemaskband()
  • GDALARGetNextUpdatedRegion() -> FUNCTION gdalargetnextupdatedregion()
  • GDALARLockBuffer() -> FUNCTION gdalarlockbuffer()
  • GDALARUnlockBuffer() -> SUBROUTINE gdalarunlockbuffer()
  • GDALSwapWords() -> SUBROUTINE gdalswapwords()
  • GDALCopyWords() -> SUBROUTINE gdalcopywords()
  • GDALCopyBits() -> SUBROUTINE gdalcopybits()
  • GDALLoadWorldFile() -> FUNCTION gdalloadworldfile()
  • GDALReadWorldFile() -> FUNCTION gdalreadworldfile()
  • GDALWriteWorldFile() -> FUNCTION gdalwriteworldfile()
  • GDALLoadOziMapFile() -> FUNCTION gdalloadozimapfile()
  • GDALReadOziMapFile() -> FUNCTION gdalreadozimapfile()
  • GDALLoadTabFile() -> FUNCTION gdalloadtabfile()
  • GDALReadTabFile() -> FUNCTION gdalreadtabfile()
  • GDALLoadRPBFile() -> FUNCTION gdalloadrpbfile()
  • GDALLoadRPCFile() -> FUNCTION gdalloadrpcfile()
  • GDALWriteRPBFile() -> FUNCTION gdalwriterpbfile()
  • GDALLoadIMDFile() -> FUNCTION gdalloadimdfile()
  • GDALWriteIMDFile() -> FUNCTION gdalwriteimdfile()
  • GDALDecToDMS() -> FUNCTION gdaldectodms()
  • GDALPackedDMSToDec() -> FUNCTION gdalpackeddmstodec()
  • GDALDecToPackedDMS() -> FUNCTION gdaldectopackeddms()
  • GDALExtractRPCInfo() -> FUNCTION gdalextractrpcinfo()
  • GDALCreateColorTable() -> FUNCTION gdalcreatecolortable()
  • GDALDestroyColorTable() -> SUBROUTINE gdaldestroycolortable()
  • GDALCloneColorTable() -> FUNCTION gdalclonecolortable()
  • GDALGetPaletteInterpretation() -> FUNCTION gdalgetpaletteinterpretation()
  • GDALGetColorEntryCount() -> FUNCTION gdalgetcolorentrycount()
  • GDALGetColorEntry() -> FUNCTION gdalgetcolorentry()
  • GDALGetColorEntryAsRGB() -> FUNCTION gdalgetcolorentryasrgb()
  • GDALSetColorEntry() -> SUBROUTINE gdalsetcolorentry()
  • GDALCreateColorRamp() -> SUBROUTINE gdalcreatecolorramp()
  • GDALSetCacheMax() -> SUBROUTINE gdalsetcachemax()
  • GDALGetCacheMax() -> FUNCTION gdalgetcachemax()
  • GDALGetCacheUsed() -> FUNCTION gdalgetcacheused()
  • GDALSetCacheMax64() -> SUBROUTINE gdalsetcachemax64()
  • GDALGetCacheMax64() -> FUNCTION gdalgetcachemax64()
  • GDALGetCacheUsed64() -> FUNCTION gdalgetcacheused64()
  • GDALFlushCacheBlock() -> FUNCTION gdalflushcacheblock()

As a general guideline, note that when a char** object is encountered in the C interface, it should usually be interfaced in Fortran by means of the fortranc::c_ptr_ptr derived type.

Other Fortran-style subroutines, functions and procedure interfaces are documented explicitely here.

For an example of application of the gdal module, please refer to the following test program, which creates a very simple gdal raster dataset, exports it on a GEOTiff file and successively reads it:

PROGRAM gdal_test
use,INTRINSIC :: iso_c_binding
USE gdal
IMPLICIT none
TYPE(gdaldriverh) :: driver
TYPE(gdaldataseth) :: ds
TYPE(gdalrasterbandh) :: band
CHARACTER(len=512) :: file
REAL(kind=c_double) :: x1, y1, x2, y2, gt(6)
INTEGER(kind=c_int) :: i1, j1, k1, i2, j2, k2, i, j, k, ierr
REAL,ALLOCATABLE :: z(:,:), z3(:,:,:), zr3(:,:,:)
! register all available gdal drivers
CALL gdalallregister()
! file name to work on
file = 'gdal_test.tiff'
! ==== How to create a gdal dataset and export it to a file ====
! get a specific driver (see e.g. gdalinfo --formats)
print*,'Getting GeoTIFF driver'
driver = gdalgetdriverbyname('GTiff'//char(0))
IF (.NOT.gdalassociated(driver)) THEN
print*,'Error getting GeoTIFF driver from gdal'
stop 1
ENDIF
! create the dataset with i1xj1 points and 1 raster band of byte data
print*,'Creating a GeoTIFF gdal dataset'
i1 = 120
j1 = 80
k1 = 3
! pass a couple of options (BIGTIFF and COMPRESS) to the create function,
! the strings must have the same legth (possibly compiler-dependent),
! they will be trimmed by the c_ptr_ptr_new function
ds = gdalcreate(driver, trim(file)//char(0), i1, j1, k1, gdt_byte, &
c_ptr_ptr_getobject(c_ptr_ptr_new((/'BIGTIFF=YES ','COMPRESS=DEFLATE'/))))
! (/('',i=1,0)/) is a trick to define a zero-length array on the fly,
! if we do not want to pass any specific option
!ds = gdalcreate(driver, TRIM(file)//CHAR(0), i1, j1, k1, gdt_byte, &
! c_ptr_ptr_getobject(c_ptr_ptr_new((/('',i=1,0)/))))
IF (.NOT.gdalassociated(ds)) THEN
print*,'Error creating a GeoTIFF dataset on file ',trim(file)
stop 1
ENDIF
! this seems to be useless here, depends on file type
print*,'Setting color interpretation to RGB'
ierr = gdalsetrastercolorinterpretation(gdalgetrasterband(ds, 1), gci_redband)
ierr = gdalsetrastercolorinterpretation(gdalgetrasterband(ds, 2), gci_greenband)
ierr = gdalsetrastercolorinterpretation(gdalgetrasterband(ds, 3), gci_blueband)
! create a dummy array of real data, we stick to integer <= 255 in
! order to fit in a byte and not to lose precision in read test
ALLOCATE(z3(i1,j1,k1))
DO k = 1, k1
DO j = 1, j1
DO i = 1, i1
z3(i,j,k) = i/2+j*mod(k,2)+(j1-j)*(1-mod(k,2))
ENDDO
ENDDO
ENDDO
! write data to dataset with the simplified Fortran interface, one
! raster band is written starting from the upper(?) left corner, the
! conversion from the type of z3 (real) to the gdt_byte type specified
! with gdalcreate is done internally by gdal
print*,'Writing data to dataset'
print*,'using the simplified Fortran interface'
ierr = gdaldatasetrasterio_f(ds, gf_write, 0, 0, z3)
IF (ierr /= 0) THEN
print*,'Error writing data to GeoTIFF dataset on file ',trim(file)
stop 1
ENDIF
CALL gdalclose(ds)
! ==== How to read a gdal dataset from a file, simplified Fortran interface ====
print*,'Opening a GeoTIFF gdal dataset for reading'
print*,'using the simplified Fortran interface'
ds = gdalopen(trim(file)//char(0), ga_readonly)
IF (.NOT.gdalassociated(ds)) THEN
print*,'Error opening dataset on file ',trim(file)
stop 1
ENDIF
print*,'Getting the affine transformation'
ierr = gdalgetgeotransform(ds, gt)
! an error is acceptable since no transformation had been defined
!IF (ierr /= 0) THEN
! PRINT*,'Error getting the affine transformation from dataset on file ',TRIM(file)
! STOP 1
!ENDIF
print*,'The affine transformation matrix is: ',gt
print*,'Getting the dataset size'
i2 = gdalgetrasterxsize(ds)
j2 = gdalgetrasterysize(ds)
k2 = gdalgetrastercount(ds)
IF (i1 /= i2 .OR. j1 /= j2) THEN
print*,'Error, wrong dataset x or y size ',i1,i2,j1,j2
stop 1
ENDIF
print*,'The x/y size of the raster is: ',i2,j2
IF (k1 /= k2) THEN
print*,'Error, wrong raster band number ',k1,k2
stop 1
ENDIF
print*,'The number of raster bands is: ',k2
! apply the affine transformation to some coordinates
! original gdal version
CALL gdalapplygeotransform(gt, 0.5_c_double, 0.5_c_double, x1, y1)
! Fortran elemental version
CALL gdalapplygeotransform_f(gt, 0.5_c_double, 0.5_c_double, x2, y2)
IF (x1 /= x2 .OR. y1 /= y2) THEN ! this check should be relaxed
print*,'Error gdal and Fortran version of GDALApplyGeoTransform'
print*,'give different results'
stop 1
ENDIF
CALL gdalapplygeotransform_f(gt, 0.5_c_double, 0.5_c_double, x1, y1)
REAL(i1, kind=c_double)-0.5_c_double, &
REAL(j1, kind=c_double)-0.5_c_double, x2, y2)
print*,'The raster coordinates of the corners are: ',x1,y1,x2,y2
ALLOCATE(z(i2,j2))
! get the first raster band
print*,'Getting a raster band'
band = gdalgetrasterband(ds, 1)
IF (.NOT.gdalassociated(band)) THEN
print*,'Error getting raster band from GeoTIFF dataset on file ',trim(file)
stop 1
ENDIF
! read data from first raster band with the simplified Fortran interface,
! raster band is read starting from the upper(?) left corner
print*,'Reading data from dataset'
ierr = gdalrasterio_f(band, gf_read, 0, 0, z)
IF (ierr /= 0) THEN
print*,'Error reading data from GeoTIFF dataset on file ',trim(file)
print*,'with simplified Fortran interface'
stop 1
ENDIF
print*,'The sum of the buffer read is: ',sum(z)
print*,'Average error after write/read process: ',&
sum(abs(z(:,:) - z3(:,:,1)))/(i2*j2)
CALL gdalclose(ds)
! ==== How to read a gdal dataset from a file, even more simplified Fortran interface ====
print*,'Opening a GeoTIFF gdal dataset for reading'
print*,'using the even more simplified Fortran interface'
ds = gdalopen(trim(file)//char(0), ga_readonly)
IF (.NOT.gdalassociated(ds)) THEN
print*,'Error opening dataset on file ',trim(file)
stop 1
ENDIF
! read data from first raster band with the even more simplified Fortran interface,
! raster band is read starting from the upper(?) left corner
print*,'Reading data from dataset'
CALL gdaldatasetsimpleread_f(ds, 5._c_double, 5._c_double, 20._c_double, 15._c_double, &
zr3, x1, y1, x2, y2)
IF (.NOT.ALLOCATED(zr3)) THEN
print*,'Error reading data from GeoTIFF dataset on file ',trim(file)
print*,'with even more simplified Fortran interface'
stop 1
ENDIF
print*,'The shape of the buffer read is: ',shape(zr3)
print*,'The sum of the buffer read is: ',sum(zr3)
print*,'The raster coordinates of the corners are: ',x1,y1,x2,y2
CALL gdalclose(ds)
DEALLOCATE(z3,z,zr3)
END PROGRAM gdal_test
Constructor for a c_ptr_ptr object.
Definition fortranc.F90:185
Interface to a Fortran version of gdalapplygeotransform working on scalars, 1-d, 2-d and 3-d arrays.
Definition gdal.F90:314
Simplified Fortran generic interface to the gdaldatasetrasterio C function.
Definition gdal.F90:343
Simplified Fortran generic interface to the gdalrasterio C function.
Definition gdal.F90:372
Utility module for supporting Fortran 2003 C language interface module.
Definition fortranc.F90:103
Fortran 2003 interface to the gdal http://www.gdal.org/ library.
Definition gdal.F90:201