main.c 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625
  1. /***************************************************************************
  2. *
  3. * MODULE: r.proj
  4. *
  5. * AUTHOR(S): Martin Schroeder
  6. * University of Heidelberg
  7. * Dept. of Geography
  8. * emes@geo0.geog.uni-heidelberg.de
  9. *
  10. * (With the help of a lot of existing GRASS sources, in
  11. * particular v.proj)
  12. *
  13. * PURPOSE: r.proj converts a map to a new geographic projection. It reads a
  14. * map from a different location, projects it and write it out
  15. * to the current location. The projected data is resampled with
  16. * one of three different methods: nearest neighbor, bilinear and
  17. * cubic convolution.
  18. *
  19. * COPYRIGHT: (C) 2001, 2011 by the GRASS Development Team
  20. *
  21. * This program is free software under the GNU General Public
  22. * License (>=v2). Read the file COPYING that comes with GRASS
  23. * for details.
  24. *
  25. * Changes
  26. * Morten Hulden <morten@untamo.net>, Aug 2000:
  27. * - aborts if input map is outside current location.
  28. * - can handle projections (conic, azimuthal etc) where
  29. * part of the map may fall into areas where south is
  30. * upward and east is leftward.
  31. * - avoids passing location edge coordinates to PROJ
  32. * (they may be invalid in some projections).
  33. * - output map will be clipped to borders of the current region.
  34. * - output map cell edges and centers will coinside with those
  35. * of the current region.
  36. * - output map resolution (unless changed explicitly) will
  37. * match (exactly) the resolution of the current region.
  38. * - if the input map is smaller than the current region, the
  39. * output map will only cover the overlapping area.
  40. * - if the input map is larger than the current region, only the
  41. * needed amount of memory will be allocated for the projection
  42. *
  43. * Bugfixes 20050328: added floor() before (int) typecasts to in avoid
  44. * asymmetrical rounding errors. Added missing offset outcellhd.ew_res/2
  45. * to initial xcoord for each row in main projection loop (we want to project
  46. * center of cell, not border).
  47. *
  48. * Glynn Clements 2006: Use G_interp_* functions, modified
  49. * version of r.proj which uses a tile cache instead of loading the
  50. * entire map into memory.
  51. * Markus Metz 2010: lanczos and lanczos fallback interpolation methods
  52. *****************************************************************************/
  53. #include <stdio.h>
  54. #include <stdlib.h>
  55. #include <string.h>
  56. #include <unistd.h>
  57. #include <grass/gis.h>
  58. #include <grass/raster.h>
  59. #include <grass/gprojects.h>
  60. #include <grass/glocale.h>
  61. #include "r.proj.h"
  62. /* modify this table to add new methods */
  63. struct menu menu[] = {
  64. {p_nearest, "nearest", "nearest neighbor"},
  65. {p_bilinear, "bilinear", "bilinear interpolation"},
  66. {p_cubic, "bicubic", "bicubic interpolation"},
  67. {p_lanczos, "lanczos", "lanczos filter"},
  68. {p_bilinear_f, "bilinear_f", "bilinear interpolation with fallback"},
  69. {p_cubic_f, "bicubic_f", "bicubic interpolation with fallback"},
  70. {p_lanczos_f, "lanczos_f", "lanczos filter with fallback"},
  71. {NULL, NULL, NULL}
  72. };
  73. static char *make_ipol_list(void);
  74. static char *make_ipol_desc(void);
  75. int main(int argc, char **argv)
  76. {
  77. char *mapname, /* ptr to name of output layer */
  78. *setname, /* ptr to name of input mapset */
  79. *ipolname; /* name of interpolation method */
  80. int fdi, /* input map file descriptor */
  81. fdo, /* output map file descriptor */
  82. method, /* position of method in table */
  83. permissions, /* mapset permissions */
  84. cell_type, /* output celltype */
  85. cell_size, /* size of a cell in bytes */
  86. row, col, /* counters */
  87. irows, icols, /* original rows, cols */
  88. orows, ocols, have_colors, /* Input map has a colour table */
  89. overwrite, /* Overwrite */
  90. curr_proj; /* output projection (see gis.h) */
  91. void *obuffer; /* buffer that holds one output row */
  92. struct cache *ibuffer; /* buffer that holds the input map */
  93. func interpolate; /* interpolation routine */
  94. double xcoord2, /* temporary x coordinates */
  95. ycoord2, /* temporary y coordinates */
  96. onorth, osouth, /* save original border coords */
  97. oeast, owest, inorth, isouth, ieast, iwest;
  98. char north_str[30], south_str[30], east_str[30], west_str[30];
  99. struct Colors colr; /* Input map colour table */
  100. struct History history;
  101. struct pj_info iproj, /* input map proj parameters */
  102. oproj, /* output map proj parameters */
  103. tproj; /* transformation parameters */
  104. struct Key_Value *in_proj_info, /* projection information of */
  105. *in_unit_info, /* input and output mapsets */
  106. *out_proj_info, *out_unit_info;
  107. struct GModule *module;
  108. struct Flag *list, /* list files in source location */
  109. *nocrop, /* don't crop output map */
  110. *print_bounds, /* print output bounds and exit */
  111. *gprint_bounds; /* same but print shell style */
  112. struct Option *imapset, /* name of input mapset */
  113. *inmap, /* name of input layer */
  114. *inlocation, /* name of input location */
  115. *outmap, /* name of output layer */
  116. *indbase, /* name of input database */
  117. *interpol, /* interpolation method */
  118. *memory, /* amount of memory for cache */
  119. *res; /* resolution of target map */
  120. #ifdef HAVE_PROJ_H
  121. struct Option *pipeline; /* name of custom PROJ pipeline */
  122. #endif
  123. struct Cell_head incellhd, /* cell header of input map */
  124. outcellhd; /* and output map */
  125. G_gisinit(argv[0]);
  126. module = G_define_module();
  127. G_add_keyword(_("raster"));
  128. G_add_keyword(_("projection"));
  129. G_add_keyword(_("transformation"));
  130. G_add_keyword(_("import"));
  131. module->description =
  132. _("Re-projects a raster map from given location to the current location.");
  133. inlocation = G_define_standard_option(G_OPT_M_LOCATION);
  134. inlocation->required = YES;
  135. inlocation->label = _("Location containing input raster map");
  136. inlocation->guisection = _("Source");
  137. imapset = G_define_standard_option(G_OPT_M_MAPSET);
  138. imapset->label = _("Mapset containing input raster map");
  139. imapset->description = _("Default: name of current mapset");
  140. imapset->guisection = _("Source");
  141. inmap = G_define_standard_option(G_OPT_R_INPUT);
  142. inmap->description = _("Name of input raster map to re-project");
  143. inmap->required = NO;
  144. inmap->guisection = _("Source");
  145. indbase = G_define_standard_option(G_OPT_M_DBASE);
  146. indbase->label = _("Path to GRASS database of input location");
  147. outmap = G_define_standard_option(G_OPT_R_OUTPUT);
  148. outmap->required = NO;
  149. outmap->description = _("Name for output raster map (default: same as 'input')");
  150. outmap->guisection = _("Target");
  151. ipolname = make_ipol_list();
  152. interpol = G_define_option();
  153. interpol->key = "method";
  154. interpol->type = TYPE_STRING;
  155. interpol->required = NO;
  156. interpol->answer = "nearest";
  157. interpol->options = ipolname;
  158. interpol->description = _("Interpolation method to use");
  159. interpol->guisection = _("Target");
  160. interpol->descriptions = make_ipol_desc();
  161. memory = G_define_option();
  162. memory->key = "memory";
  163. memory->type = TYPE_INTEGER;
  164. memory->required = NO;
  165. memory->answer = "300";
  166. memory->label = _("Maximum memory to be used (in MB)");
  167. memory->description = _("Cache size for raster rows");
  168. res = G_define_option();
  169. res->key = "resolution";
  170. res->type = TYPE_DOUBLE;
  171. res->required = NO;
  172. res->description = _("Resolution of output raster map");
  173. res->guisection = _("Target");
  174. #ifdef HAVE_PROJ_H
  175. pipeline = G_define_option();
  176. pipeline->key = "pipeline";
  177. pipeline->type = TYPE_STRING;
  178. pipeline->required = NO;
  179. pipeline->description = _("PROJ pipeline for coordinate transformation");
  180. #endif
  181. list = G_define_flag();
  182. list->key = 'l';
  183. list->description = _("List raster maps in input mapset and exit");
  184. list->guisection = _("Print");
  185. nocrop = G_define_flag();
  186. nocrop->key = 'n';
  187. nocrop->description = _("Do not perform region cropping optimization");
  188. print_bounds = G_define_flag();
  189. print_bounds->key = 'p';
  190. print_bounds->description =
  191. _("Print input map's bounds in the current projection and exit");
  192. print_bounds->guisection = _("Print");
  193. gprint_bounds = G_define_flag();
  194. gprint_bounds->key = 'g';
  195. gprint_bounds->description =
  196. _("Print input map's bounds in the current projection and exit (shell style)");
  197. gprint_bounds->guisection = _("Print");
  198. /* The parser checks if the map already exists in current mapset,
  199. we switch out the check and do it
  200. in the module after the parser */
  201. overwrite = G_check_overwrite(argc, argv);
  202. if (G_parser(argc, argv))
  203. exit(EXIT_FAILURE);
  204. /* get the method */
  205. for (method = 0; (ipolname = menu[method].name); method++)
  206. if (strcmp(ipolname, interpol->answer) == 0)
  207. break;
  208. if (!ipolname)
  209. G_fatal_error(_("<%s=%s> unknown %s"),
  210. interpol->key, interpol->answer, interpol->key);
  211. interpolate = menu[method].method;
  212. mapname = outmap->answer ? outmap->answer : inmap->answer;
  213. if (mapname && !list->answer && !overwrite &&
  214. !print_bounds->answer && !gprint_bounds->answer &&
  215. G_find_raster(mapname, G_mapset()))
  216. G_fatal_error(_("option <%s>: <%s> exists. To overwrite, use the --overwrite flag"),
  217. "output", mapname);
  218. setname = imapset->answer ? imapset->answer : G_store(G_mapset());
  219. if (strcmp(inlocation->answer, G_location()) == 0 &&
  220. (!indbase->answer || strcmp(indbase->answer, G_gisdbase()) == 0))
  221. #if 0
  222. G_fatal_error(_("Input and output locations can not be the same"));
  223. #else
  224. G_warning(_("Input and output locations are the same"));
  225. #endif
  226. G_get_window(&outcellhd);
  227. if (gprint_bounds->answer && !print_bounds->answer)
  228. print_bounds->answer = gprint_bounds->answer;
  229. curr_proj = G_projection();
  230. /* Get projection info for output mapset */
  231. if ((out_proj_info = G_get_projinfo()) == NULL)
  232. G_fatal_error(_("Unable to get projection info of output raster map"));
  233. if ((out_unit_info = G_get_projunits()) == NULL)
  234. G_fatal_error(_("Unable to get projection units of output raster map"));
  235. if (pj_get_kv(&oproj, out_proj_info, out_unit_info) < 0)
  236. G_fatal_error(_("Unable to get projection key values of output raster map"));
  237. /* Change the location */
  238. G_create_alt_env();
  239. G_setenv_nogisrc("GISDBASE", indbase->answer ? indbase->answer : G_gisdbase());
  240. G_setenv_nogisrc("LOCATION_NAME", inlocation->answer);
  241. permissions = G_mapset_permissions(setname);
  242. if (permissions < 0) /* can't access mapset */
  243. G_fatal_error(_("Mapset <%s> in input location <%s> - %s"),
  244. setname, inlocation->answer,
  245. permissions == 0 ? _("permission denied")
  246. : _("not found"));
  247. /* if requested, list the raster maps in source location - MN 5/2001 */
  248. if (list->answer) {
  249. int i;
  250. char **srclist;
  251. G_verbose_message(_("Checking location <%s> mapset <%s>"),
  252. inlocation->answer, setname);
  253. srclist = G_list(G_ELEMENT_RASTER, G_getenv_nofatal("GISDBASE"),
  254. G_getenv_nofatal("LOCATION_NAME"), setname);
  255. for (i = 0; srclist[i]; i++) {
  256. fprintf(stdout, "%s\n", srclist [i]);
  257. }
  258. fflush(stdout);
  259. exit(EXIT_SUCCESS); /* leave r.proj after listing */
  260. }
  261. if (!inmap->answer)
  262. G_fatal_error(_("Required parameter <%s> not set"), inmap->key);
  263. if (!G_find_raster(inmap->answer, setname))
  264. G_fatal_error(_("Raster map <%s> in location <%s> in mapset <%s> not found"),
  265. inmap->answer, inlocation->answer, setname);
  266. /* Read input map colour table */
  267. have_colors = Rast_read_colors(inmap->answer, setname, &colr);
  268. /* Get projection info for input mapset */
  269. if ((in_proj_info = G_get_projinfo()) == NULL)
  270. G_fatal_error(_("Unable to get projection info of input map"));
  271. /* apparently the +over switch must be set in the input projection,
  272. * not the output latlon projection */
  273. if (curr_proj == PROJECTION_LL)
  274. G_set_key_value("+over", "defined", in_proj_info);
  275. if ((in_unit_info = G_get_projunits()) == NULL)
  276. G_fatal_error(_("Unable to get projection units of input map"));
  277. if (pj_get_kv(&iproj, in_proj_info, in_unit_info) < 0)
  278. G_fatal_error(_("Unable to get projection key values of input map"));
  279. tproj.def = NULL;
  280. #ifdef HAVE_PROJ_H
  281. if (pipeline->answer) {
  282. tproj.def = G_store(pipeline->answer);
  283. }
  284. #endif
  285. if (GPJ_init_transform(&iproj, &oproj, &tproj) < 0)
  286. G_fatal_error(_("Unable to initialize coordinate transformation"));
  287. G_free_key_value(in_proj_info);
  288. G_free_key_value(in_unit_info);
  289. G_free_key_value(out_proj_info);
  290. G_free_key_value(out_unit_info);
  291. if (G_verbose() > G_verbose_std())
  292. pj_print_proj_params(&iproj, &oproj);
  293. /* this call causes r.proj to read the entire map into memeory */
  294. Rast_get_cellhd(inmap->answer, setname, &incellhd);
  295. Rast_set_input_window(&incellhd);
  296. if (G_projection() == PROJECTION_XY)
  297. G_fatal_error(_("Unable to work with unprojected data (xy location)"));
  298. /* Save default borders so we can show them later */
  299. inorth = incellhd.north;
  300. isouth = incellhd.south;
  301. ieast = incellhd.east;
  302. iwest = incellhd.west;
  303. irows = incellhd.rows;
  304. icols = incellhd.cols;
  305. onorth = outcellhd.north;
  306. osouth = outcellhd.south;
  307. oeast = outcellhd.east;
  308. owest = outcellhd.west;
  309. orows = outcellhd.rows;
  310. ocols = outcellhd.cols;
  311. if (print_bounds->answer) {
  312. G_message(_("Input map <%s@%s> in location <%s>:"),
  313. inmap->answer, setname, inlocation->answer);
  314. outcellhd.north = -1e9;
  315. outcellhd.south = 1e9;
  316. outcellhd.east = -1e9;
  317. outcellhd.west = 1e9;
  318. bordwalk_edge(&incellhd, &outcellhd, &iproj, &oproj, &tproj, PJ_FWD);
  319. inorth = outcellhd.north;
  320. isouth = outcellhd.south;
  321. ieast = outcellhd.east;
  322. iwest = outcellhd.west;
  323. G_format_northing(inorth, north_str, curr_proj);
  324. G_format_northing(isouth, south_str, curr_proj);
  325. G_format_easting(ieast, east_str, curr_proj);
  326. G_format_easting(iwest, west_str, curr_proj);
  327. if (gprint_bounds->answer) {
  328. fprintf(stdout, "n=%s s=%s w=%s e=%s rows=%d cols=%d\n",
  329. north_str, south_str, west_str, east_str, irows, icols);
  330. }
  331. else {
  332. fprintf(stdout, "Source cols: %d\n", icols);
  333. fprintf(stdout, "Source rows: %d\n", irows);
  334. fprintf(stdout, "Local north: %s\n", north_str);
  335. fprintf(stdout, "Local south: %s\n", south_str);
  336. fprintf(stdout, "Local west: %s\n", west_str);
  337. fprintf(stdout, "Local east: %s\n", east_str);
  338. }
  339. exit(EXIT_SUCCESS);
  340. }
  341. /* Cut non-overlapping parts of input map */
  342. if (!nocrop->answer)
  343. bordwalk(&outcellhd, &incellhd, &iproj, &oproj, &tproj, PJ_INV);
  344. /* Add 2 cells on each side for bilinear/cubic & future interpolation methods */
  345. /* (should probably be a factor based on input and output resolution) */
  346. incellhd.north += 2 * incellhd.ns_res;
  347. incellhd.east += 2 * incellhd.ew_res;
  348. incellhd.south -= 2 * incellhd.ns_res;
  349. incellhd.west -= 2 * incellhd.ew_res;
  350. if (incellhd.north > inorth)
  351. incellhd.north = inorth;
  352. if (incellhd.east > ieast)
  353. incellhd.east = ieast;
  354. if (incellhd.south < isouth)
  355. incellhd.south = isouth;
  356. if (incellhd.west < iwest)
  357. incellhd.west = iwest;
  358. Rast_set_input_window(&incellhd);
  359. /* And switch back to original location */
  360. G_switch_env();
  361. /* Adjust borders of output map */
  362. if (!nocrop->answer)
  363. bordwalk(&incellhd, &outcellhd, &iproj, &oproj, &tproj, PJ_FWD);
  364. #if 0
  365. outcellhd.west = outcellhd.south = HUGE_VAL;
  366. outcellhd.east = outcellhd.north = -HUGE_VAL;
  367. for (row = 0; row < incellhd.rows; row++) {
  368. ycoord1 = Rast_row_to_northing((double)(row + 0.5), &incellhd);
  369. for (col = 0; col < incellhd.cols; col++) {
  370. xcoord1 = Rast_col_to_easting((double)(col + 0.5), &incellhd);
  371. if (GPJ_transform(&iproj, &oproj, &tproj, PJ_FWD,
  372. &xcoord1, &ycoord1, NULL) < 0)
  373. G_fatal_error(_("Error in %s"), "GPJ_transform()");
  374. if (xcoord1 > outcellhd.east)
  375. outcellhd.east = xcoord1;
  376. if (ycoord1 > outcellhd.north)
  377. outcellhd.north = ycoord1;
  378. if (xcoord1 < outcellhd.west)
  379. outcellhd.west = xcoord1;
  380. if (ycoord1 < outcellhd.south)
  381. outcellhd.south = ycoord1;
  382. }
  383. }
  384. #endif
  385. if (res->answer != NULL) /* set user defined resolution */
  386. outcellhd.ns_res = outcellhd.ew_res = atof(res->answer);
  387. G_adjust_Cell_head(&outcellhd, 0, 0);
  388. Rast_set_output_window(&outcellhd);
  389. G_message(" ");
  390. G_message(_("Input:"));
  391. G_message(_("Cols: %d (%d)"), incellhd.cols, icols);
  392. G_message(_("Rows: %d (%d)"), incellhd.rows, irows);
  393. G_message(_("North: %f (%f)"), incellhd.north, inorth);
  394. G_message(_("South: %f (%f)"), incellhd.south, isouth);
  395. G_message(_("West: %f (%f)"), incellhd.west, iwest);
  396. G_message(_("East: %f (%f)"), incellhd.east, ieast);
  397. G_message(_("EW-res: %f"), incellhd.ew_res);
  398. G_message(_("NS-res: %f"), incellhd.ns_res);
  399. G_message(" ");
  400. G_message(_("Output:"));
  401. G_message(_("Cols: %d (%d)"), outcellhd.cols, ocols);
  402. G_message(_("Rows: %d (%d)"), outcellhd.rows, orows);
  403. G_message(_("North: %f (%f)"), outcellhd.north, onorth);
  404. G_message(_("South: %f (%f)"), outcellhd.south, osouth);
  405. G_message(_("West: %f (%f)"), outcellhd.west, owest);
  406. G_message(_("East: %f (%f)"), outcellhd.east, oeast);
  407. G_message(_("EW-res: %f"), outcellhd.ew_res);
  408. G_message(_("NS-res: %f"), outcellhd.ns_res);
  409. G_message(" ");
  410. /* open and read the relevant parts of the input map and close it */
  411. G_switch_env();
  412. Rast_set_input_window(&incellhd);
  413. fdi = Rast_open_old(inmap->answer, setname);
  414. cell_type = Rast_get_map_type(fdi);
  415. ibuffer = readcell(fdi, memory->answer);
  416. Rast_close(fdi);
  417. G_switch_env();
  418. Rast_set_output_window(&outcellhd);
  419. if (strcmp(interpol->answer, "nearest") == 0) {
  420. fdo = Rast_open_new(mapname, cell_type);
  421. obuffer = (CELL *) Rast_allocate_output_buf(cell_type);
  422. }
  423. else {
  424. fdo = Rast_open_fp_new(mapname);
  425. cell_type = FCELL_TYPE;
  426. obuffer = (FCELL *) Rast_allocate_output_buf(cell_type);
  427. }
  428. cell_size = Rast_cell_size(cell_type);
  429. xcoord2 = outcellhd.west + (outcellhd.ew_res / 2);
  430. ycoord2 = outcellhd.north - (outcellhd.ns_res / 2);
  431. G_important_message(_("Projecting..."));
  432. for (row = 0; row < outcellhd.rows; row++) {
  433. /* obufptr = obuffer */;
  434. G_percent(row, outcellhd.rows - 1, 2);
  435. #if 0
  436. /* parallelization does not always work,
  437. * segfaults in the interpolation functions
  438. * can happen */
  439. #pragma omp parallel for schedule (static)
  440. #endif
  441. for (col = 0; col < outcellhd.cols; col++) {
  442. void *obufptr = (void *)((const unsigned char *)obuffer + col * cell_size);
  443. double xcoord1 = xcoord2 + (col) * outcellhd.ew_res;
  444. double ycoord1 = ycoord2;
  445. /* project coordinates in output matrix to */
  446. /* coordinates in input matrix */
  447. if (GPJ_transform(&iproj, &oproj, &tproj, PJ_INV,
  448. &xcoord1, &ycoord1, NULL) < 0) {
  449. G_warning(_("Error in %s"), "GPJ_transform()");
  450. Rast_set_null_value(obufptr, 1, cell_type);
  451. }
  452. else {
  453. /* convert to row/column indices of input matrix */
  454. /* column index in input matrix */
  455. double col_idx = (xcoord1 - incellhd.west) / incellhd.ew_res;
  456. /* row index in input matrix */
  457. double row_idx = (incellhd.north - ycoord1) / incellhd.ns_res;
  458. /* and resample data point */
  459. interpolate(ibuffer, obufptr, cell_type,
  460. col_idx, row_idx, &incellhd);
  461. }
  462. /* obufptr = G_incr_void_ptr(obufptr, cell_size); */
  463. }
  464. Rast_put_row(fdo, obuffer, cell_type);
  465. xcoord2 = outcellhd.west + (outcellhd.ew_res / 2);
  466. ycoord2 -= outcellhd.ns_res;
  467. }
  468. Rast_close(fdo);
  469. release_cache(ibuffer);
  470. if (have_colors > 0) {
  471. Rast_write_colors(mapname, G_mapset(), &colr);
  472. Rast_free_colors(&colr);
  473. }
  474. Rast_short_history(mapname, "raster", &history);
  475. Rast_command_history(&history);
  476. Rast_write_history(mapname, &history);
  477. G_done_msg(" ");
  478. exit(EXIT_SUCCESS);
  479. }
  480. char *make_ipol_list(void)
  481. {
  482. int size = 0;
  483. int i;
  484. char *buf;
  485. for (i = 0; menu[i].name; i++)
  486. size += strlen(menu[i].name) + 1;
  487. buf = G_malloc(size);
  488. *buf = '\0';
  489. for (i = 0; menu[i].name; i++) {
  490. if (i)
  491. strcat(buf, ",");
  492. strcat(buf, menu[i].name);
  493. }
  494. return buf;
  495. }
  496. char *make_ipol_desc(void)
  497. {
  498. int size = 0;
  499. int i;
  500. char *buf;
  501. for (i = 0; menu[i].name; i++)
  502. size += strlen(menu[i].name) + strlen(menu[i].text) + 2;
  503. buf = G_malloc(size);
  504. *buf = '\0';
  505. for (i = 0; menu[i].name; i++) {
  506. if (i)
  507. strcat(buf, ";");
  508. strcat(buf, menu[i].name);
  509. strcat(buf, ";");
  510. strcat(buf, menu[i].text);
  511. }
  512. return buf;
  513. }