GRASS GIS 7 Programmer's Manual  7.9.dev(2021)-e5379bbd7
gdal.c
Go to the documentation of this file.
1 /*!
2  \file lib/raster/gdal.c
3 
4  \brief Raster Library - Utilization of GDAL library.
5 
6  (C) 2010 by the GRASS Development Team
7 
8  This program is free software under the GNU General Public License
9  (>=v2). Read the file COPYING that comes with GRASS for details.
10 
11  \author Glynn Clements
12 */
13 
14 #include <stdlib.h>
15 #include <string.h>
16 #include <unistd.h>
17 
18 #include <grass/config.h>
19 #include <grass/gis.h>
20 #include <grass/raster.h>
21 #include <grass/gprojects.h>
22 #include <grass/glocale.h>
23 
24 #include "R.h"
25 
26 #ifndef HAVE_GDAL
27 #undef GDAL_LINK
28 #endif
29 
30 #ifdef GDAL_LINK
31 
32 #ifdef GDAL_DYNAMIC
33 # if defined(__unix) || defined(__unix__)
34 # include <dlfcn.h>
35 # endif
36 # ifdef _WIN32
37 # include <windows.h>
38 # endif
39 #endif
40 
41 static void CPL_STDCALL(*pGDALAllRegister) (void);
42 static void CPL_STDCALL(*pGDALClose) (GDALDatasetH);
43 static GDALRasterBandH CPL_STDCALL(*pGDALGetRasterBand) (GDALDatasetH, int);
44 static GDALDatasetH CPL_STDCALL(*pGDALOpen) (const char *pszFilename,
45  GDALAccess eAccess);
46 static CPLErr CPL_STDCALL(*pGDALRasterIO) (GDALRasterBandH hRBand,
47  GDALRWFlag eRWFlag, int nDSXOff,
48  int nDSYOff, int nDSXSize,
49  int nDSYSize, void *pBuffer,
50  int nBXSize, int nBYSize,
51  GDALDataType eBDataType,
52  int nPixelSpace, int nLineSpace);
53 static GDALDriverH CPL_STDCALL(*pGDALGetDriverByName) (const char *);
54 static const char *CPL_STDCALL(*pGDALGetMetadataItem) (GDALMajorObjectH,
55  const char *,
56  const char *);
57 static GDALDatasetH CPL_STDCALL(*pGDALCreate) (GDALDriverH hDriver,
58  const char *, int, int, int,
59  GDALDataType, char **);
60 static GDALDatasetH CPL_STDCALL(*pGDALCreateCopy) (GDALDriverH, const char *,
61  GDALDatasetH, int, char **,
62  GDALProgressFunc, void *);
63 static CPLErr CPL_STDCALL(*pGDALSetRasterNoDataValue) (GDALRasterBandH,
64  double);
65 static CPLErr CPL_STDCALL(*pGDALSetGeoTransform) (GDALDatasetH, double *);
66 static CPLErr CPL_STDCALL(*pGDALSetProjection) (GDALDatasetH, const char *);
67 static const char *CPL_STDCALL(*pGDALGetDriverShortName) (GDALDriverH);
68 static GDALDriverH CPL_STDCALL(*pGDALGetDatasetDriver) (GDALDatasetH);
69 
70 #if GDAL_DYNAMIC
71 # if defined(__unix) && !defined(__unix__)
72 # define __unix__ __unix
73 # endif
74 
75 static void *library_h;
76 
77 static void *get_symbol(const char *name)
78 {
79  void *sym;
80 
81 # ifdef __unix__
82  sym = dlsym(library_h, name);
83 # endif
84 # ifdef _WIN32
85  sym = GetProcAddress((HINSTANCE) library_h, name);
86 # endif
87 
88  if (!sym)
89  G_fatal_error(_("Unable to locate symbol <%s>"), name);
90 
91  return sym;
92 }
93 
94 static void try_load_library(const char *name)
95 {
96 # ifdef __unix__
97  library_h = dlopen(name, RTLD_NOW);
98 # endif
99 # ifdef _WIN32
100  library_h = LoadLibrary(name);
101 # endif
102 }
103 
104 static void load_library(void)
105 {
106  static const char *const candidates[] = {
107 # ifdef __unix__
108  "libgdal.so.26", /* GDAL 3.0 */
109  "libgdal.so.20",
110  "libgdal.so.1",
111  "libgdal.1.1.so",
112  "gdal.1.0.so",
113  "gdal.so.1.0",
114  "libgdal.so",
115  "libgdal1.6.0.so",
116  "libgdal1.7.0.so",
117 # endif
118 # ifdef _WIN32
119  "gdal300.dll",
120  "gdal204.dll",
121  "gdal203.dll",
122  "gdal202.dll",
123  "gdal201.dll",
124  "gdal200.dll",
125  "gdal111.dll",
126  "gdal110.dll",
127  "gdal19.dll",
128  "gdal18.dll",
129  "gdal17.dll",
130  "gdal16.dll",
131  "gdal15.dll",
132  "gdal11.dll",
133  "gdal.1.0.dll",
134  "libgdal-1.dll",
135  "gdal.dll",
136 # endif
137  NULL
138  };
139  int i;
140 
141  for (i = 0; candidates[i]; i++) {
142  try_load_library(candidates[i]);
143  if (library_h) {
144  G_debug(3, "found %s", candidates[i]);
145  return;
146  }
147  }
148 
149  G_fatal_error(_("Unable to load GDAL library"));
150 }
151 
152 static void init_gdal(void)
153 {
154  load_library();
155 
156 # if defined(_WIN32) && !defined(_WIN64)
157  pGDALAllRegister = get_symbol("_GDALAllRegister@0");
158  pGDALOpen = get_symbol("_GDALOpen@8");
159  pGDALClose = get_symbol("_GDALClose@4");
160  pGDALGetRasterBand = get_symbol("_GDALGetRasterBand@8");
161  pGDALRasterIO = get_symbol("_GDALRasterIO@48");
162  pGDALGetDriverByName = get_symbol("_GDALGetDriverByName@4");
163  pGDALGetMetadataItem = get_symbol("_GDALGetMetadataItem@12");
164  pGDALCreate = get_symbol("_GDALCreate@28");
165  pGDALCreateCopy = get_symbol("_GDALCreateCopy@28");
166  pGDALSetRasterNoDataValue = get_symbol("_GDALSetRasterNoDataValue@12");
167  pGDALSetGeoTransform = get_symbol("_GDALSetGeoTransform@8");
168  pGDALSetProjection = get_symbol("_GDALSetProjection@8");
169  pGDALGetDriverShortName = get_symbol("_GDALGetDriverShortName@4");
170  pGDALGetDatasetDriver = get_symbol("_GDALGetDatasetDriver@4");
171 #else
172  pGDALAllRegister = get_symbol("GDALAllRegister");
173  pGDALOpen = get_symbol("GDALOpen");
174  pGDALClose = get_symbol("GDALClose");
175  pGDALGetRasterBand = get_symbol("GDALGetRasterBand");
176  pGDALRasterIO = get_symbol("GDALRasterIO");
177  pGDALGetDriverByName = get_symbol("GDALGetDriverByName");
178  pGDALGetMetadataItem = get_symbol("GDALGetMetadataItem");
179  pGDALCreate = get_symbol("GDALCreate");
180  pGDALCreateCopy = get_symbol("GDALCreateCopy");
181  pGDALSetRasterNoDataValue = get_symbol("GDALSetRasterNoDataValue");
182  pGDALSetGeoTransform = get_symbol("GDALSetGeoTransform");
183  pGDALSetProjection = get_symbol("GDALSetProjection");
184  pGDALGetDriverShortName = get_symbol("GDALGetDriverShortName");
185  pGDALGetDatasetDriver = get_symbol("GDALGetDatasetDriver");
186 #endif
187 }
188 
189 #else /* GDAL_DYNAMIC */
190 
191 static void init_gdal(void)
192 {
193  pGDALAllRegister = &GDALAllRegister;
194  pGDALOpen = &GDALOpen;
195  pGDALClose = &GDALClose;
196  pGDALGetRasterBand = &GDALGetRasterBand;
197  pGDALRasterIO = &GDALRasterIO;
198  pGDALGetDriverByName = &GDALGetDriverByName;
199  pGDALGetMetadataItem = &GDALGetMetadataItem;
200  pGDALCreate = &GDALCreate;
201  pGDALCreateCopy = &GDALCreateCopy;
202  pGDALSetRasterNoDataValue = &GDALSetRasterNoDataValue;
203  pGDALSetGeoTransform = &GDALSetGeoTransform;
204  pGDALSetProjection = &GDALSetProjection;
205  pGDALGetDriverShortName = &GDALGetDriverShortName;
206  pGDALGetDatasetDriver = &GDALGetDatasetDriver;
207 }
208 
209 #endif /* GDAL_DYNAMIC */
210 
211 #endif /* GDAL_LINK */
212 
213 /*!
214  \brief Initialization
215 
216  Register all GDAL drivers.
217 */
218 void Rast_init_gdal(void)
219 {
220 #ifdef GDAL_LINK
221  static int initialized;
222 
223  if (G_is_initialized(&initialized))
224  return;
225 
226  init_gdal();
227  (*pGDALAllRegister) ();
228  G_initialize_done(&initialized);
229 #endif
230 }
231 
232 /*!
233  \brief Get GDAL link settings for given raster map
234 
235  \param name map name
236  \param mapset name of mapset
237 
238  \return pointer to GDAL_link structure
239  \return NULL if link not found
240 */
241 struct GDAL_link *Rast_get_gdal_link(const char *name, const char *mapset)
242 {
243 #ifdef GDAL_LINK
244  GDALDatasetH data;
245  GDALRasterBandH band;
246  GDALDataType type;
247  RASTER_MAP_TYPE req_type;
248 #endif
249  const char *filename;
250  int band_num;
251  struct GDAL_link *gdal;
252  RASTER_MAP_TYPE map_type;
253  FILE *fp;
254  struct Key_Value *key_val;
255  const char *p;
256  DCELL null_val;
257  int hflip, vflip;
258 
259  if (!G_find_raster2(name, mapset))
260  return NULL;
261 
262  map_type = Rast_map_type(name, mapset);
263  if (map_type < 0)
264  return NULL;
265 
266  fp = G_fopen_old_misc("cell_misc", "gdal", name, mapset);
267  if (!fp)
268  return NULL;
269  key_val = G_fread_key_value(fp);
270  fclose(fp);
271 
272  if (!key_val)
273  return NULL;
274 
275  filename = G_find_key_value("file", key_val);
276  if (!filename)
277  return NULL;
278 
279  p = G_find_key_value("band", key_val);
280  if (!p)
281  return NULL;
282  band_num = atoi(p);
283  if (!band_num)
284  return NULL;
285 
286  p = G_find_key_value("null", key_val);
287  if (!p)
288  return NULL;
289  /* atof on windows can not read "nan" and returns 0 instead */
290  if (strcmp(p, "none") == 0 ||
291  G_strcasecmp(p, "nan") == 0 || G_strcasecmp(p, "-nan") == 0) {
292  Rast_set_d_null_value(&null_val, 1);
293  }
294  else
295  null_val = atof(p);
296 
297  hflip = G_find_key_value("hflip", key_val) ? 1 : 0;
298  vflip = G_find_key_value("vflip", key_val) ? 1 : 0;
299 
300 #ifdef GDAL_LINK
301  p = G_find_key_value("type", key_val);
302  if (!p)
303  return NULL;
304  type = atoi(p);
305 
306  switch (type) {
307  case GDT_Byte:
308  case GDT_Int16:
309  case GDT_UInt16:
310  case GDT_Int32:
311  case GDT_UInt32:
312  req_type = CELL_TYPE;
313  break;
314  case GDT_Float32:
315  req_type = FCELL_TYPE;
316  break;
317  case GDT_Float64:
318  req_type = DCELL_TYPE;
319  break;
320  default:
321  return NULL;
322  }
323 
324  if (req_type != map_type)
325  return NULL;
326 
327  Rast_init_gdal();
328 
329  data = (*pGDALOpen) (filename, GA_ReadOnly);
330  if (!data)
331  return NULL;
332 
333  band = (*pGDALGetRasterBand) (data, band_num);
334  if (!band) {
335  (*pGDALClose) (data);
336  return NULL;
337  }
338 #endif
339 
340  gdal = G_calloc(1, sizeof(struct GDAL_link));
341 
342  gdal->filename = G_store(filename);
343  gdal->band_num = band_num;
344  gdal->null_val = null_val;
345  gdal->hflip = hflip;
346  gdal->vflip = vflip;
347 #ifdef GDAL_LINK
348  gdal->data = data;
349  gdal->band = band;
350  gdal->type = type;
351 #endif
352 
353  return gdal;
354 }
355 
356 struct GDAL_Options
357 {
358  const char *dir;
359  const char *ext;
360  const char *format;
361  char **options;
362 };
363 
364 static struct state
365 {
366  int initialized;
367  struct GDAL_Options opts;
368  struct Key_Value *projinfo, *projunits, *projepsg;
369  char *srswkt;
370 } state;
371 
372 static struct state *st = &state;
373 
374 static void read_gdal_options(void)
375 {
376  FILE *fp;
377  struct Key_Value *key_val;
378  const char *p;
379 
380  fp = G_fopen_old("", "GDAL", G_mapset());
381  if (!fp)
382  G_fatal_error(_("Unable to open GDAL file"));
383  key_val = G_fread_key_value(fp);
384  fclose(fp);
385 
386  p = G_find_key_value("directory", key_val);
387  if (!p)
388  p = "gdal";
389  if (*p == '/') {
390  st->opts.dir = G_store(p);
391  }
392  else {
393  char path[GPATH_MAX];
394 
395  G_file_name(path, p, "", G_mapset());
396  st->opts.dir = G_store(path);
397  if (access(path, 0) != 0)
399  }
400 
401  p = G_find_key_value("extension", key_val);
402  st->opts.ext = G_store(p ? p : "");
403 
404  p = G_find_key_value("format", key_val);
405  st->opts.format = G_store(p ? p : "GTiff");
406 
407  p = G_find_key_value("options", key_val);
408  st->opts.options = p ? G_tokenize(p, ",") : NULL;
409 
410  G_free_key_value(key_val);
411 }
412 
413 /*!
414  \brief Create GDAL settings for given raster map
415 
416  \param name map name
417  \param map_type map type (CELL, FCELL, DCELL)
418 
419  \return pointer to allocated GDAL_link structure
420  \return NULL on error
421 */
422 struct GDAL_link *Rast_create_gdal_link(const char *name,
423  RASTER_MAP_TYPE map_type)
424 {
425 #ifdef GDAL_LINK
426  char path[GPATH_MAX];
427  GDALDriverH driver;
428  double transform[6];
429  struct GDAL_link *gdal;
430  FILE *fp;
431  struct Key_Value *key_val;
432  char buf[32];
433 
435 
436  Rast_init_gdal();
437 
438  if (!G_is_initialized(&st->initialized)) {
439  read_gdal_options();
440  st->projinfo = G_get_projinfo();
441  st->projunits = G_get_projunits();
442  st->projepsg = G_get_projepsg();
443  if (st->projinfo && st->projunits)
444  st->srswkt = GPJ_grass_to_wkt2(st->projinfo, st->projunits,
445  st->projepsg, 0, 0);
446  G_initialize_done(&st->initialized);
447  }
448 
449  gdal = G_calloc(1, sizeof(struct GDAL_link));
450 
451  sprintf(path, "%s/%s%s", st->opts.dir, name, st->opts.ext);
452  gdal->filename = G_store(path);
453  gdal->band_num = 1;
454  gdal->hflip = 0;
455  gdal->vflip = 0;
456 
457  switch (map_type) {
458  case CELL_TYPE:
459  switch (R__.nbytes) {
460  case 1:
461  gdal->type = GDT_Byte;
462  gdal->null_val = (DCELL) 0xFF;
463  break;
464  case 2:
465  gdal->type = GDT_UInt16;
466  gdal->null_val = (DCELL) 0xFFFF;
467  break;
468  case 3:
469  case 4:
470  gdal->type = GDT_Int32;
471  gdal->null_val = (DCELL) 0x80000000U;
472  break;
473  }
474  break;
475  case FCELL_TYPE:
476  gdal->type = GDT_Float32;
477  Rast_set_d_null_value(&gdal->null_val, 1);
478  break;
479  case DCELL_TYPE:
480  gdal->type = GDT_Float64;
481  Rast_set_d_null_value(&gdal->null_val, 1);
482  break;
483  default:
484  G_fatal_error(_("Invalid map type <%d>"), map_type);
485  break;
486  }
487 
488  driver = (*pGDALGetDriverByName) (st->opts.format);
489  if (!driver)
490  G_fatal_error(_("Unable to get <%s> driver"), st->opts.format);
491 
492  /* Does driver support GDALCreate ? */
493  if ((*pGDALGetMetadataItem) (driver, GDAL_DCAP_CREATE, NULL)) {
494  gdal->data =
495  (*pGDALCreate)(driver, gdal->filename,
497  1, gdal->type, st->opts.options);
498  if (!gdal->data)
499  G_fatal_error(_("Unable to create <%s> dataset using <%s> driver"),
500  name, st->opts.format);
501  }
502  /* If not - create MEM driver for intermediate dataset.
503  * Check if raster can be created at all (with GDALCreateCopy) */
504  else if ((*pGDALGetMetadataItem) (driver, GDAL_DCAP_CREATECOPY, NULL)) {
505  GDALDriverH mem_driver;
506 
507  G_message(_("Driver <%s> does not support direct writing. "
508  "Using MEM driver for intermediate dataset."),
509  st->opts.format);
510 
511  mem_driver = (*pGDALGetDriverByName) ("MEM");
512  if (!mem_driver)
513  G_fatal_error(_("Unable to get in-memory raster driver"));
514 
515  gdal->data =
516  (*pGDALCreate)(mem_driver, "",
518  1, gdal->type, st->opts.options);
519  if (!gdal->data)
520  G_fatal_error(_("Unable to create <%s> dataset using memory driver"),
521  name);
522  }
523  else
524  G_fatal_error(_("Driver <%s> does not support creating rasters"),
525  st->opts.format);
526 
527  gdal->band = (*pGDALGetRasterBand) (gdal->data, gdal->band_num);
528 
529  (*pGDALSetRasterNoDataValue) (gdal->band, gdal->null_val);
530 
531  /* Set Geo Transform */
532  transform[0] = R__.wr_window.west;
533  transform[1] = R__.wr_window.ew_res;
534  transform[2] = 0.0;
535  transform[3] = R__.wr_window.north;
536  transform[4] = 0.0;
537  transform[5] = -R__.wr_window.ns_res;
538 
539  if ((*pGDALSetGeoTransform) (gdal->data, transform) >= CE_Failure)
540  G_warning(_("Unable to set geo transform"));
541 
542  if (st->srswkt)
543  if ((*pGDALSetProjection) (gdal->data, st->srswkt) == CE_Failure)
544  G_warning(_("Unable to set projection"));
545 
546  fp = G_fopen_new_misc("cell_misc", "gdal", name);
547  if (!fp)
548  G_fatal_error(_("Unable to create cell_misc/%s/gdal file"), name);
549 
550  key_val = G_create_key_value();
551 
552  G_set_key_value("file", gdal->filename, key_val);
553 
554  sprintf(buf, "%d", gdal->band_num);
555  G_set_key_value("band", buf, key_val);
556 
557  sprintf(buf, "%.22g", gdal->null_val);
558  G_set_key_value("null", buf, key_val);
559 
560  sprintf(buf, "%d", gdal->type);
561  G_set_key_value("type", buf, key_val);
562 
563  if (G_fwrite_key_value(fp, key_val) < 0)
564  G_fatal_error(_("Error writing cell_misc/%s/gdal file"), name);
565 
566  G_free_key_value(key_val);
567 
568  fclose(fp);
569 
570  return gdal;
571 #else
572  return NULL;
573 #endif
574 }
575 
576 /*!
577  \brief Close existing GDAL link
578 
579  \param gdal pointer to GDAL_link to be closed
580 */
581 void Rast_close_gdal_link(struct GDAL_link *gdal)
582 {
583 #ifdef GDAL_LINK
584  (*pGDALClose) (gdal->data);
585 #endif
586  G_free(gdal->filename);
587  G_free(gdal);
588 }
589 
590 /*!
591  \brief Close existing GDAL link and write out data
592 
593  \param gdal pointer to GDAL_link to be closed
594 
595  \return 1 on success
596  \return -1 on failure
597 */
599 {
600  int stat = 1;
601 
602 #ifdef GDAL_LINK
603  GDALDriverH src_drv = (*pGDALGetDatasetDriver) (gdal->data);
604 
605  if (G_strcasecmp((*pGDALGetDriverShortName) (src_drv), "MEM") == 0) {
606  GDALDriverH dst_drv = (*pGDALGetDriverByName) (st->opts.format);
607  GDALDatasetH dst =
608  (*pGDALCreateCopy) (dst_drv, gdal->filename, gdal->data, FALSE,
609  st->opts.options, NULL, NULL);
610 
611  if (!dst) {
612  G_warning(_("Unable to create output file <%s> using driver <%s>"),
613  gdal->filename, st->opts.format);
614  stat = -1;
615  }
616  (*pGDALClose) (dst);
617  }
618 
619  (*pGDALClose) (gdal->data);
620 
621 #endif
622  G_free(gdal->filename);
623  G_free(gdal);
624 
625  return stat;
626 }
627 
628 #ifdef GDAL_LINK
629 /*!
630  \brief Input/output function for GDAL links
631 
632  See GDAL's RasterIO for details.
633 */
634 CPLErr Rast_gdal_raster_IO(GDALRasterBandH band, GDALRWFlag rw_flag,
635  int x_off, int y_off, int x_size, int y_size,
636  void *buffer, int buf_x_size, int buf_y_size,
637  GDALDataType buf_type, int pixel_size,
638  int line_size)
639 {
640  return (*pGDALRasterIO) (band, rw_flag, x_off, y_off, x_size, y_size,
641  buffer, buf_x_size, buf_y_size, buf_type,
642  pixel_size, line_size);
643 }
644 #endif
#define CELL_TYPE
Definition: raster.h:11
char * G_file_name(char *, const char *, const char *, const char *)
Builds full path names to GIS data files.
Definition: file_name.c:38
void void void void G_fatal_error(const char *,...) __attribute__((format(printf
struct Key_Value * G_get_projinfo(void)
Gets projection information for location.
Definition: get_projinfo.c:61
Definition: R.h:88
CPLErr Rast_gdal_raster_IO(GDALRasterBandH, GDALRWFlag, int, int, int, int, void *, int, int, GDALDataType, int, int)
int Rast_close_gdal_write_link(struct GDAL_link *gdal)
Close existing GDAL link and write out data.
Definition: gdal.c:598
double west
Extent coordinates (west)
Definition: gis.h:464
double DCELL
Definition: gis.h:603
struct Key_Value * G_fread_key_value(FILE *)
Read key/values pairs from file.
Definition: key_value2.c:49
void G_free(void *)
Free allocated memory.
Definition: gis/alloc.c:149
char * dst
Definition: lz4.h:599
char ** G_tokenize(const char *, const char *)
Tokenize string.
Definition: gis/token.c:48
#define NULL
Definition: ccmath.h:32
void Rast__init_window(void)
void G_initialize_done(int *)
Definition: counter.c:76
#define G_calloc(m, n)
Definition: defs/gis.h:113
struct Cell_head wr_window
Definition: R.h:100
struct Key_Value * G_get_projunits(void)
Gets units information for location.
Definition: get_projinfo.c:32
void G_message(const char *,...) __attribute__((format(printf
struct state * st
Definition: parser.c:104
#define DCELL_TYPE
Definition: raster.h:13
double north
Extent coordinates (north)
Definition: gis.h:458
struct GDAL_link * Rast_create_gdal_link(const char *name, RASTER_MAP_TYPE map_type)
Create GDAL settings for given raster map.
Definition: gdal.c:422
int int G_strcasecmp(const char *, const char *)
String compare ignoring case (upper or lower)
Definition: strings.c:47
#define FALSE
Definition: gis.h:63
char * GPJ_grass_to_wkt2(const struct Key_Value *, const struct Key_Value *, const struct Key_Value *, int, int)
Converts a GRASS co-ordinate system representation to WKT style. EPSG code is preferred if available...
Definition: convert.c:152
FILE * G_fopen_old_misc(const char *, const char *, const char *, const char *)
open a database file for reading
Definition: open_misc.c:210
const struct driver * driver
Definition: driver/init.c:25
void G_set_key_value(const char *, const char *, struct Key_Value *)
Set value for given key.
Definition: key_value1.c:38
#define GPATH_MAX
Definition: gis.h:170
FILE * G_fopen_old(const char *, const char *, const char *)
Open a database file for reading.
Definition: gis/open.c:253
Definition: gis.h:501
int G_fwrite_key_value(FILE *, const struct Key_Value *)
Write key/value pairs to file.
Definition: key_value2.c:25
int cols
Number of columns for 2D data.
Definition: gis.h:431
int G_is_initialized(int *)
Definition: counter.c:59
FILE * G_fopen_new_misc(const char *, const char *, const char *)
open a new database file
Definition: open_misc.c:182
double ns_res
Resolution - north to south cell size for 2D data.
Definition: gis.h:452
const char * G_mapset(void)
Get current mapset name.
Definition: gis/mapset.c:33
void G_warning(const char *,...) __attribute__((format(printf
void Rast_init_gdal(void)
Initialization.
Definition: gdal.c:218
Definition: path.h:16
#define _(str)
Definition: glocale.h:10
int RASTER_MAP_TYPE
Definition: raster.h:25
int G_make_mapset_element(const char *)
Create element in the current mapset.
Definition: mapset_msc.c:38
#define FCELL_TYPE
Definition: raster.h:12
const char * G_find_raster2(const char *, const char *)
Find a raster map (look but don&#39;t touch)
Definition: find_rast.c:76
char * G_store(const char *)
Copy string to allocated memory.
Definition: strings.c:87
void Rast_close_gdal_link(struct GDAL_link *gdal)
Close existing GDAL link.
Definition: gdal.c:581
struct Key_Value * G_create_key_value(void)
Allocate and initialize Key_Value structure.
Definition: key_value1.c:23
int nbytes
Definition: R.h:94
const char * name
Definition: named_colr.c:7
double ew_res
Resolution - east to west cell size for 2D data.
Definition: gis.h:448
int rows
Number of rows for 2D data.
Definition: gis.h:427
struct state state
Definition: parser.c:103
int G_debug(int, const char *,...) __attribute__((format(printf
void Rast_set_d_null_value(DCELL *, int)
To set a number of DCELL raster values to NULL.
Definition: null_val.c:155
void G_free_key_value(struct Key_Value *)
Free allocated Key_Value structure.
Definition: key_value1.c:103
struct Key_Value * G_get_projepsg(void)
Gets EPSG information for the current location.
Definition: get_projinfo.c:102
struct GDAL_link * Rast_get_gdal_link(const char *name, const char *mapset)
Get GDAL link settings for given raster map.
Definition: gdal.c:241
const char * G_find_key_value(const char *, const struct Key_Value *)
Find given key (case sensitive)
Definition: key_value1.c:84
RASTER_MAP_TYPE Rast_map_type(const char *, const char *)
Determine raster data type.
Definition: raster/open.c:880